Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JPulseFinder_utils.hh File Reference
#include <string>
#include <iostream>
#include <sstream>
#include <iomanip>
#include "JDAQ/JDAQ.hh"
#include "JDAQ/JDAQTimeslice.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JHitToolkit.hh"
#include "JTrigger/JTriggerParameters.hh"
#include "JTrigger/JTimesliceRouter.hh"
#include "TH1D.h"
#include "TH2D.h"

Go to the source code of this file.

Functions

int get_offset (JPMT pmt, JModule NB, double c_w, IO options)
 Calculates an offset for a combinataion nanobeacon-pmt that allows to define the histogram range in a time window containing the corresponding pulse. More...
 
vector< vector< TH2D * > > allocate_histograms (JDetector detector, IO options, double c_w, vector< int > top_pmts, vector< int > bottom_pmts, int max_distance)
 Allocates histograms to store the tot vs time hit distribution for all the combinations nanobeacon - pmt. More...
 
void fill_histograms (vector< vector< TH2D * > > &histograms, IO options, JDetector detector, vector< int > top_pmts, vector< int > bottom_pmts, int max_distance)
 Fills the hit time vs ToT histograms for all the combinations nanobeacon - pmt. More...
 
vector< vector< TH1D * > > project_ToT (vector< vector< TH2D * > > ToT_vs_time)
 Projects the hit time vs ToT histograms for all the combinations nanobeacon - pmt on the Y axis. More...
 
vector< vector< TH1D * > > project_time (vector< vector< TH2D * > > ToT_vs_time)
 Projects the hit time vs ToT histograms for all the combinations nanobeacon - pmt on the X axis. More...
 

Function Documentation

int get_offset ( JPMT  pmt,
JModule  NB,
double  c_w,
IO  options 
)
inline

Calculates an offset for a combinataion nanobeacon-pmt that allows to define the histogram range in a time window containing the corresponding pulse.

Author
rgruiz
Parameters
pmtThe target pmt
NBThe source module
c_wThe speed light in water
optionsuser options containing the stagger and pulse delay information
Returns
The corresponding offset to apply in the definition of the hit time histogram range.

Definition at line 53 of file JPulseFinder_utils.hh.

53  {
54 
55  double distance_pmt_nb = (pmt.getPosition() - NB.getPMT(0).getPosition()).getLength() ;
56 
57  double time_of_flight = distance_pmt_nb / c_w ;
58 
59  double delta_t0 = pmt.getT0() - NB.getPMT(0).getT0() ;
60 
61  double delta_traw = time_of_flight - delta_t0 ;
62 
63  double stagger_delay = 16 * options.stagger_16ns * (NB.getFloor() - 1) ;
64 
65  double overall_delay = 16 * options.overall_delay_16ns ;
66 
67  int total_offset = (int)round( delta_traw + stagger_delay + overall_delay);
68 
69  return total_offset ;
70 
71 }
int stagger_16ns
int overall_delay_16ns
int getFloor() const
Get floor number.
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:129
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:141
double getT0() const
Get time offset.
vector< vector < TH2D* > > allocate_histograms ( JDetector  detector,
IO  options,
double  c_w,
vector< int >  top_pmts,
vector< int >  bottom_pmts,
int  max_distance 
)
inline

Allocates histograms to store the tot vs time hit distribution for all the combinations nanobeacon - pmt.

Parameters
detectorThe detector containing the DU to be analyzed
optionsuser options containing the string number and the list of used nanobeacons
c_wThe speed light in water
top_pmtsList of FPGA channels for pmts in the upper hemisphere
bottom_pmtsList of FPGA channels for pmts in the lower hemisphere.
max_distanceMaximum distance between modules
Returns
A table of TH2D histograms defined in the correct hit time ranges.

Definition at line 86 of file JPulseFinder_utils.hh.

86  {
87 
88  int nFloors = getNumberOfFloors(detector) ;
89 
90  int nHistosPerFloor = top_pmts.size() + max_distance * bottom_pmts.size() ;
91 
92  vector < vector < TH2D* > > histograms (nFloors , vector<TH2D*> (nHistosPerFloor, NULL)) ;
93 
94  char name[250] ;
95 
96  for (int tgt_floor = 1 ; tgt_floor <= nFloors ; ++tgt_floor){
97 
98  JModule tgt_module = detector.getModule( JModuleLocation(options.string_number , tgt_floor) ) ;
99 
100  int index = 0 ;
101 
102  for(auto & top_pmt : top_pmts){
103 
104  int total_offset = get_offset(tgt_module.getPMT(top_pmt) , tgt_module , c_w , options) ;
105 
106  sprintf(name , "S%iF%iPMT%i_beaconS%iF%i" , tgt_module.getString() , tgt_module.getFloor() , top_pmt , tgt_module.getString() , tgt_module.getFloor());
107 
108  histograms[tgt_floor - 1][index] = new TH2D( name , NULL , 400 , total_offset - 0.5 , 400 + total_offset - 0.5 , 256 , -0.5 , 255.5 ) ;
109 
110  histograms[tgt_floor - 1][index]->SetOption("colz") ;
111 
112  index += 1 ;
113  }
114 
115 
116  for (int src_floor = tgt_floor - max_distance ; src_floor < tgt_floor ; ++src_floor){
117 
118  if(src_floor > 0){
119 
120  JModule src_module = detector.getModule( JModuleLocation (options.string_number , src_floor) ) ;
121 
122  for (auto & bottom_pmt : bottom_pmts){
123 
124  int total_offset = get_offset(tgt_module.getPMT(bottom_pmt) , src_module , c_w , options) ;
125 
126  sprintf(name , "S%iF%iPMT%i_beaconS%iF%i" , tgt_module.getString() , tgt_module.getFloor() , bottom_pmt , src_module.getString() , src_module.getFloor());
127 
128  histograms[tgt_floor - 1][index] = new TH2D( name , NULL , 400 , total_offset - 0.5 , 400 + total_offset - 0.5 , 256 , -0.5 , 255.5 ) ;
129 
130  histograms[tgt_floor - 1][index]->SetOption("colz") ;
131 
132  index += 1 ;
133 
134  }
135 
136  }
137 
138  }
139 
140  }
141 
142  return histograms ;
143 
144 }
int string_number
Data structure for a composite optical module.
Definition: JModule.hh:47
int getFloor() const
Get floor number.
int getString() const
Get string number.
int getNumberOfFloors(const JDetector &detector)
Get number of floors.
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:141
const JModule & getModule(const JModuleAddress &address) const
Get module parameters.
Definition: JDetector.hh:180
int get_offset(JPMT pmt, JModule NB, double c_w, IO options)
Calculates an offset for a combinataion nanobeacon-pmt that allows to define the histogram range in a...
Logical location of module.
void fill_histograms ( vector< vector< TH2D * > > &  histograms,
IO  options,
JDetector  detector,
vector< int >  top_pmts,
vector< int >  bottom_pmts,
int  max_distance 
)
inline

Fills the hit time vs ToT histograms for all the combinations nanobeacon - pmt.

Parameters
histogramsThe vector of histograms to be filled
optionsUser options containing the string number and the list files to be read
detectorThe detector containing the DU to be analyzed
top_pmtsList of FPGA channels for pmts in the upper hemisphere
bottom_pmtsList of FPGA channels for pmts in the lower hemisphere.
max_distanceMaximum distance between modules

Definition at line 159 of file JPulseFinder_utils.hh.

159  {
160 
161  JTriggerParameters parameters ;
162 
163  unsigned int pulse_period = 16 * options.pulse_period_16ns ;
164 
166 
167  JModuleRouter moduleRouter(detector) ;
168 
169  int slicecounter = 0 ;
170 
171  while( inputFile.hasNext() ) {
172 
173  slicecounter ++ ;
174 
175  JDAQTimesliceL0 slice = *(inputFile.next()) ;
176 
177  JTimesliceRouter timesliceRouter(parameters.numberOfBins) ;
178 
179  timesliceRouter.configure(slice);
180 
181  for (int i = 0 ; i < getNumberOfFloors(detector) ; i++){
182 
183  const JModule& module = detector.getModule(JModuleLocation(options.string_number , i+1)) ;
184 
185  if(timesliceRouter.hasSuperFrame(module.getID())){
186 
187  JDAQSuperFrame SuperFrame = timesliceRouter.getSuperFrame(module.getID()) ;
188 
189  vector<JDAQHit> buffer[getNumberOfPMTs(module)];
190 
191  for (JDAQSuperFrame::const_iterator hit = SuperFrame.begin(); hit != SuperFrame.end(); ++hit) {
192 
193  buffer[hit->getPMT()].push_back(*hit) ;
194 
195  }
196 
197  int i = 0 ;
198 
199  for (auto & top_pmt : top_pmts){
200 
201  for (int hit = 0 ; hit < (int)buffer[top_pmt].size() - 1 ; hit++){
202 
203  double tot = buffer[top_pmt][hit].getToT() ;
204 
205  double raw_time = buffer[top_pmt][hit].getT() ;
206 
207  //double t0 = module.getPMT(buffer[top_pmt][hit].getPMT()).getT0() ;
208 
209  //double calibrated_time = raw_time + t0 ;
210 
211  double mod_double = raw_time - pulse_period * int(raw_time/pulse_period) ;
212 
213  histograms[module.getFloor()-1][i]->Fill(mod_double , tot) ;
214 
215  }
216 
217  i++ ;
218 
219  }
220 
221 
222  i = 0 ;
223 
224  for (auto & bottom_pmt : bottom_pmts){
225 
226  for (int hit = 0 ; hit < (int)buffer[bottom_pmt].size() - 1 ; hit++){
227 
228  double tot = buffer[bottom_pmt][hit].getToT() ;
229 
230  double raw_time = buffer[bottom_pmt][hit].getT() ;
231 
232  //double t0 = module.getPMT(buffer[bottom_pmt][hit].getPMT()).getT0() ;
233 
234  //double calibrated_time = raw_time + t0 ;
235 
236  double mod_double = raw_time - pulse_period * int(raw_time/pulse_period) ;
237 
238  for (int j = 1 ; j <= max_distance ; j++){
239 
240  if (module.getFloor() - j <=0) continue ;
241 
242  histograms[module.getFloor() - 1][top_pmts.size() + (j-1)*bottom_pmts.size() + i]->Fill(mod_double , tot) ;
243 
244  }
245 
246  }
247 
248  i++ ;
249 
250  }
251 
252  } // end of loop over JDAQTimeSlices
253 
254  }
255 
256  }
257 
258 }
unsigned int pulse_period_16ns
Data structure for all trigger parameters.
int string_number
Data structure for a composite optical module.
Definition: JModule.hh:47
Router for direct addressing of module data in detector data structure.
int numberOfBins
number of bins for lookup table of timeslice
int getNumberOfPMTs(const JModule &module)
Get number of PMTs.
JMultipleFileScanner ifnames
int getFloor() const
Get floor number.
Hit data structure.
Definition: JDAQHit.hh:36
Router for fast addressing of hits in KM3NETDAQ::JDAQTimeslice data structure as a function of the op...
const_iterator begin() const
Definition: JDAQFrame.hh:139
int getID() const
Get identifier.
Definition: JObjectID.hh:54
int getNumberOfFloors(const JDetector &detector)
Get number of floors.
const JModule & getModule(const JModuleAddress &address) const
Get module parameters.
Definition: JDetector.hh:180
General purpose class for object reading from a list of file names.
Data frame of one optical module.
Timeslice data structure for L0 data.
Logical location of module.
const_iterator end() const
Definition: JDAQFrame.hh:140
vector< vector < TH1D* > > project_ToT ( vector< vector< TH2D * > >  ToT_vs_time)

Projects the hit time vs ToT histograms for all the combinations nanobeacon - pmt on the Y axis.

Parameters
ToT_vs_timeThe vector of histograms to be projected
Returns
A vector with ToT distribiutions.

Definition at line 270 of file JPulseFinder_utils.hh.

270  {
271 
272  vector < vector < TH1D* > > ToT (ToT_vs_time.size() , vector < TH1D* > (ToT_vs_time[0].size() , NULL)) ;
273 
274  for (int i=0 ; i<(int)ToT_vs_time.size() ; i++){
275 
276  for (int j=0 ; j<(int)ToT_vs_time[0].size() ; j++){
277 
278  if (ToT_vs_time[i][j]==NULL) continue ;
279 
280  TH1D* p_y = ToT_vs_time[i][j]->ProjectionY((string("ToT_") + ToT_vs_time[i][j]->GetName()).c_str() , 1 , ToT_vs_time[i][j]->GetNbinsX()-1) ;
281 
282  ToT[i][j] = p_y ;
283 
284  }
285 
286  }
287 
288  return ToT ;
289 
290 }
vector< vector < TH1D* > > project_time ( vector< vector< TH2D * > >  ToT_vs_time)
inline

Projects the hit time vs ToT histograms for all the combinations nanobeacon - pmt on the X axis.

Parameters
ToT_vs_timeThe vector of histograms to be projected
Returns
A vector with hit time distribiutions.

Definition at line 300 of file JPulseFinder_utils.hh.

300  {
301 
302  vector < vector < TH1D* > > time (ToT_vs_time.size() , vector < TH1D* > (ToT_vs_time[0].size() , NULL)) ;
303 
304  for (int i=0 ; i<(int)ToT_vs_time.size() ; i++){
305 
306  for (int j=0 ; j<(int)ToT_vs_time[0].size() ; j++){
307 
308  if (ToT_vs_time[i][j]==NULL) continue ;
309 
310  TH1D* p_x = ToT_vs_time[i][j]->ProjectionX((string("htime_") + ToT_vs_time[i][j]->GetName()).c_str() , 1 , ToT_vs_time[i][j]->GetNbinsY()-1) ;
311 
312  time[i][j] = p_x ;
313 
314  }
315 
316  }
317 
318  return time ;
319 
320 }