Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_RF_online.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_RF_online.cc
4 // Created: Wed Nov 8 11:58:09 EST 2015
5 // Creator: pmatt (on Linux stan.jlab.org 2.6.32-279.11.1.el6.x86_64 x86_64)
6 
8 
9 // Routine used to create our JEventProcessor
10 extern "C"
11 {
12  void InitPlugin(JApplication *app)
13  {
14  InitJANAPlugin(app);
15  app->AddProcessor(new JEventProcessor_RF_online());
16  }
17 }
18 
19 //DOCUMENTATION: https://halldweb.jlab.org/wiki/index.php/RF_Calibration
20 
22 {
23  //This constant should be fixed for the lifetime of GlueX. If it ever changes, move it into the CCDB.
24  dRFSignalPeriod = 1000.0/499.0; //2.004008016
25  double locDeltaTRangeMax = 2.2;
26  unsigned int locNumDeltaTBins = 1100;
27 
28  dRFSignalSystems.push_back(SYS_FDC); dRFSignalSystems.push_back(SYS_PSC);
29  dRFSignalSystems.push_back(SYS_TAGH); dRFSignalSystems.push_back(SYS_TOF);
30 
35 
36  dMaxDeltaTHits[SYS_TOF] = 15;
40 
41  string locHistName, locHistTitle;
42 
43  gDirectory->Cd("/");
44  new TDirectoryFile("RF", "RF");
45  gDirectory->cd("RF");
46 
47  //roc info consistency:
48  //compare roc infos: delta-t of each to avg-t-exlcuding-itself //separate histograms for each delta
49  //will create histograms in evnt() method
50  dROCTIDirectory = new TDirectoryFile("ROCTIs", "ROCTIs");
51 
52  //# rf signals per event
53  new TDirectoryFile("NumRFSignals", "NumRFSignals");
54  gDirectory->cd("NumRFSignals");
55  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
56  {
57  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
58  string locSystemName = SystemName(locSystem);
59 
60  locHistName = string("Num") + locSystemName + string("RFSignals");
61  locHistTitle = string("RF_") + locSystemName + string(";Num RF Signals");
62  dHistMap_NumSignals[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 20, -0.5, 19.5);
63  }
64  gDirectory->cd("..");
65 
66  //rf hits found
67  new TDirectoryFile("RFHitsFound", "RFHitsFound");
68  gDirectory->cd("RFHitsFound");
69  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
70  {
71  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
72  string locSystemName = SystemName(locSystem);
73 
74  locHistName = locSystemName + string("RFHitsFound");
75  locHistTitle = string("RF_") + locSystemName + string(";Which RF Hits Were Found");
76  dHistMap_RFHitsFound[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 20, -0.5, 19.5);
77  }
78  gDirectory->cd("..");
79 
80  //num rf hits missing
81  new TDirectoryFile("NumRFHitsMissing", "NumRFHitsMissing");
82  gDirectory->cd("NumRFHitsMissing");
83  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
84  {
85  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
86  string locSystemName = SystemName(locSystem);
87 
88  locHistName = locSystemName + string("NumRFHitsMissing");
89  locHistTitle = string("RF_") + locSystemName + string(";# RF Hits Missing");
90  dHistMap_NumRFHitsMissing[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 20, -0.5, 19.5);
91  }
92  gDirectory->cd("..");
93 
94  //rf times: //delta-t to first time
95  new TDirectoryFile("DeltaT_RF_FirstTime", "DeltaT_RF_FirstTime");
96  gDirectory->cd("DeltaT_RF_FirstTime");
97  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
98  {
99  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
100  string locSystemName = SystemName(locSystem);
101 
102  locHistName = locSystemName + string("RF_FirstTimeDeltaT");
103  locHistTitle = string("RF_") + locSystemName + string("_TDC;#Deltat (Subsequent Times - First Time) (ns)");
104  dHistMap_RFFirstTimeDeltaT[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 1000000, 0.0, 5000.0);
105  }
106  gDirectory->cd("..");
107 
108  //all delta-t's:
109  new TDirectoryFile("DeltaT_All", "DeltaT_All");
110  gDirectory->cd("DeltaT_All");
111  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
112  {
113  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
114  string locSystemName = SystemName(locSystem);
115 
116  for(size_t loc_j = 0; loc_j < dMaxDeltaTHits[locSystem]; ++loc_j)
117  {
118  for(size_t loc_k = loc_j; loc_k < dMaxDeltaTHits[locSystem]; ++loc_k)
119  {
120  pair<size_t, size_t> locTimePair(loc_j, loc_k);
121  ostringstream locHistNameStream, locHistTitleStream;
122  locHistNameStream << locSystemName << "RF_DeltaT_" << loc_j << "_" << loc_k;
123  double locHistRangeMin = -200.0 + (loc_k - loc_j)*dRFSignalPeriod*dRFSamplingFactor[locSystem];
124  locHistTitleStream << "RF_" << locSystemName << "_TDC;" << "#Deltat: Hit " << loc_k << " - Hit " << loc_j << " (ns)";
125  dHistMap_AdjacentRFDeltaTs[locSystem][locTimePair] = new TH1I(locHistNameStream.str().c_str(), locHistTitleStream.str().c_str(), 80000, locHistRangeMin, locHistRangeMin + 400.0);
126  }
127  }
128  }
129  gDirectory->cd("..");
130 
131  //rf signal period
132  new TDirectoryFile("RF_SignalPeriod", "RF_SignalPeriod");
133  gDirectory->cd("RF_SignalPeriod");
134  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
135  {
136  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
137  string locSystemName = SystemName(locSystem);
138 
139  locHistName = locSystemName + string("RF_SignalPeriod");
140  locHistTitle = string("RF_") + locSystemName + string(";RF Signal Period (ns)");
141  dHistMap_RFSignalPeriod[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 1200, 2.001, 2.007);
142  }
143  gDirectory->cd("..");
144 
145  //RF beam bunch period
146  new TDirectoryFile("RF_BeamBunchPeriod", "RF_BeamBunchPeriod");
147  gDirectory->cd("RF_BeamBunchPeriod");
148  locHistName = "RFBeamBunchPeriod";
149  locHistTitle = ";TAGH #Deltat (Within Same Counter) (ns)";
150  dHist_RFBeamBunchPeriod = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 1000, 0.0, 200.0);
151  gDirectory->cd("..");
152 
153  //resolutions: compare each rf to itself
154  new TDirectoryFile("DeltaT_RF_Itself", "DeltaT_RF_Itself");
155  gDirectory->cd("DeltaT_RF_Itself");
156  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
157  {
158  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
159  string locSystemName = SystemName(locSystem);
160 
161  locHistName = locSystemName + string("RF_SelfDeltaT");
162  locHistTitle = string("RF_") + locSystemName + string(";#Deltat (First Pair) (ns)");
163  dHistMap_SelfResolution[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), locNumDeltaTBins, -1.0*locDeltaTRangeMax, locDeltaTRangeMax);
164  }
165  gDirectory->cd("..");
166 
167  //absolute resolutions: compare each rf to each rf
168  new TDirectoryFile("AbsoluteDeltaT_RF_OtherRFs", "AbsoluteDeltaT_RF_OtherRFs");
169  gDirectory->cd("AbsoluteDeltaT_RF_OtherRFs");
170  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
171  {
172  string locSystemName1 = SystemName(dRFSignalSystems[loc_i]);
173  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
174  {
175  string locSystemName2 = SystemName(dRFSignalSystems[loc_j]);
176  locHistName = string("RFDeltaT_") + locSystemName1 + string("_") + locSystemName2;
177  locHistTitle = string(";#Deltat (RF_") + locSystemName1 + string(" - RF_") + locSystemName2 + string(")");
178  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(dRFSignalSystems[loc_i], dRFSignalSystems[loc_j]);
179  dHistMap_AbsoluteRFRFDeltaTs[locSystemPair] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 4000000, -4000.0, 4000.0);
180  }
181  }
182  gDirectory->cd("..");
183 
184  //resolutions: compare each rf to each rf
185  new TDirectoryFile("DeltaT_RF_OtherRFs", "DeltaT_RF_OtherRFs");
186  gDirectory->cd("DeltaT_RF_OtherRFs");
187  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
188  {
189  string locSystemName1 = SystemName(dRFSignalSystems[loc_i]);
190  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
191  {
192  string locSystemName2 = SystemName(dRFSignalSystems[loc_j]);
193  locHistName = string("RFDeltaT_") + locSystemName1 + string("_") + locSystemName2;
194  locHistTitle = string(";#Deltat (RF_") + locSystemName1 + string(" - RF_") + locSystemName2 + string(")");
195  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(dRFSignalSystems[loc_i], dRFSignalSystems[loc_j]);
196  dHistMap_RFRFDeltaTs[locSystemPair] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), locNumDeltaTBins, -1.0*locDeltaTRangeMax, locDeltaTRangeMax);
197  }
198  }
199  gDirectory->cd("..");
200 
201  //resolutions: compare avg of each rf to avg of each rf
202  //calib time offsets
203  new TDirectoryFile("AverageDeltaT_RF_OtherRFs", "AverageDeltaT_RF_OtherRFs");
204  gDirectory->cd("AverageDeltaT_RF_OtherRFs");
205  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
206  {
207  string locSystemName1 = SystemName(dRFSignalSystems[loc_i]);
208  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
209  {
210  string locSystemName2 = SystemName(dRFSignalSystems[loc_j]);
211  locHistName = string("RFDeltaT_") + locSystemName1 + string("_") + locSystemName2;
212  locHistTitle = string(";#Deltat_{#mu} (RF_") + locSystemName1 + string(" - RF_") + locSystemName2 + string(")");
213  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(dRFSignalSystems[loc_i], dRFSignalSystems[loc_j]);
214  dHistMap_AverageRFRFDeltaTs[locSystemPair] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), locNumDeltaTBins, -1.0*locDeltaTRangeMax, locDeltaTRangeMax);
215  }
216  }
217  gDirectory->cd("..");
218 
219  //compare each rf to tagger hodoscope times
220  new TDirectoryFile("DeltaT_RF_TAGH", "DeltaT_RF_TAGH");
221  gDirectory->cd("DeltaT_RF_TAGH");
222  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
223  {
224  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
225  string locSystemName = SystemName(locSystem);
226 
227  locHistName = locSystemName + string("RF_TaggerDeltaT");
228  locHistTitle = string(";#Deltat (RF_") + locSystemName + string(" - TAGH) (ns)");
229  dHistMap_RFTaggerDeltaT[locSystem] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), locNumDeltaTBins, -1.0*locDeltaTRangeMax, locDeltaTRangeMax);
230  }
231  gDirectory->cd("..");
232 
233  gDirectory->cd(".."); //END (get back to base folder)
234 
235  return NOERROR;
236 }
237 
238 jerror_t JEventProcessor_RF_online::brun(JEventLoop* locEventLoop, int32_t runnumber)
239 {
240  // This is called whenever the run number changes
241 
242  return NOERROR;
243 }
244 
245 jerror_t JEventProcessor_RF_online::evnt(JEventLoop* locEventLoop, uint64_t eventnumber)
246 {
247  vector<const DCODAROCInfo*> locCODAROCInfos;
248  locEventLoop->Get(locCODAROCInfos);
249 
250  const DTTabUtilities* locTTabUtilities = NULL;
251  locEventLoop->GetSingle(locTTabUtilities);
252 
253  vector<const DRFTDCDigiTime*> locRFTDCDigiTimes;
254  locEventLoop->Get(locRFTDCDigiTimes);
255 
256  vector<const DTAGHHit*> locTAGHHits;
257  locEventLoop->Get(locTAGHHits);
258 
259  auto dRFTimeFactory = static_cast<DRFTime_factory*>(locEventLoop->GetFactory("DRFTime"));
260 
261  //Convert TDCs to Times
262  //Use std::set: times are NOT necessarily in order (for high-resolution mode, times interleaved between different internal channels)
263  map<DetectorSystem_t, set<double> > locRFTimes;
264  map<DetectorSystem_t, set<double> >::iterator locTDCIterator;
265  for(size_t loc_i = 0; loc_i < locRFTDCDigiTimes.size(); ++loc_i)
266  {
267  DetectorSystem_t locSystem = locRFTDCDigiTimes[loc_i]->dSystem;
268  double locRFTime = dRFTimeFactory->Convert_TDCToTime(locRFTDCDigiTimes[loc_i], locTTabUtilities);
269  locRFTimes[locSystem].insert(locRFTime);
270  }
271 
272  //Get Total ROC Info Time
273  uint64_t locTotalROCTimeStamp = 0;
274  for(size_t loc_i = 0; loc_i < locCODAROCInfos.size(); ++loc_i)
275  locTotalROCTimeStamp += locCODAROCInfos[loc_i]->timestamp;
276 
277  //Average RF Times
278  map<DetectorSystem_t, double> locAverageRFTimes;
279  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
280  {
281  set<double>& locRFTimeSet = locTDCIterator->second;
282 
283  double locFirstTime = dRFTimeFactory->Step_TimeToNearInputTime(*(locRFTimeSet.begin()), 0.0);
284  double locAverageTime = locFirstTime;
285  set<double>::iterator locSetIterator = locRFTimeSet.begin();
286  for(++locSetIterator; locSetIterator != locRFTimeSet.end(); ++locSetIterator)
287  locAverageTime += dRFTimeFactory->Step_TimeToNearInputTime(*locSetIterator, locFirstTime);
288  locAverageTime /= double(locRFTimeSet.size());
289  locAverageRFTimes[locTDCIterator->first] = locAverageTime;
290  }
291 
292  //organize TAGH hits
293  map<int, set<double> > locTAGHHitMap; //double: TAGH tdc times (set allows time sorting)
294  for(size_t loc_i = 0; loc_i < locTAGHHits.size(); ++loc_i)
295  {
296  if(locTAGHHits[loc_i]->has_TDC)
297  locTAGHHitMap[locTAGHHits[loc_i]->counter_id].insert(locTAGHHits[loc_i]->time_tdc);
298  }
299 
300  //collect TAGH delta-t's for RF beam bunch histogram
301  map<int, set<double> >::iterator locTAGHIterator = locTAGHHitMap.begin();
302  vector<double> locTAGHDeltaTs;
303  for(; locTAGHIterator != locTAGHHitMap.end(); ++locTAGHIterator)
304  {
305  set<double>& locCounterHits = locTAGHIterator->second;
306  if(locCounterHits.size() < 2)
307  continue;
308  set<double>::iterator locSetIterator = locCounterHits.begin();
309  set<double>::iterator locPreviousSetIterator = locSetIterator;
310  for(++locSetIterator; locSetIterator != locCounterHits.end(); ++locSetIterator, ++locPreviousSetIterator)
311  locTAGHDeltaTs.push_back(*locSetIterator - *locPreviousSetIterator);
312  }
313 
314  //CALC RF SIGNAL SAMPLING RATES
315  //Set by EPICS guis, but can back-compute rather than carry around more data
316  //Compute for each event, so don't need to acquire locks for accessing/modifying member variable (processor member variables are NOT thread-local)
317  //Doesn't quite work!: Hits seem to be missing sometimes
318 /*
319  map<DetectorSystem_t, double> locRFSamplingRateMap;
320  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
321  {
322  DetectorSystem_t locSystem = locTDCIterator->first;
323  set<double>& locRFTimeSet = locTDCIterator->second;
324  if(locRFTimeSet.size() < 2)
325  continue;
326 
327  double locFirstDeltaT = *(++locRFTimeSet.begin()) - *(locRFTimeSet.begin());
328  double locRFSamplingRate = locFirstDeltaT/dRFSignalPeriod; //is a number like "128.0" so is actually an inverse rate (rate = 1/128)
329 
330  //get the nearest even integer
331  int locRFSamplingRateInt = int(locRFSamplingRate); //round down
332  if(locRFSamplingRateInt % 2 == 1) //if true: rounded to odd: wrong direction
333  locRFSamplingRateInt += (locRFSamplingRate > double(locRFSamplingRateInt)) ? 1 : -1; //if rounded down/up, increase/decrease (was wrong direction)
334  //else rounded down to even number: OK (rounding up would have been to an odd number, so rounding down was closest int)
335 
336  locRFSamplingRateMap[locSystem] = 1.0/double(locRFSamplingRateInt); //convert from int -> double, then from inverse-rate to rate (i.e. is now something like 1/128.0)
337  set<double>::iterator locSetIterator = locRFTimeSet.begin();
338  for(++locSetIterator; locSetIterator != locRFTimeSet.end(); ++locSetIterator)
339  dHistMap_RFFirstTimeDeltaT[locTDCIterator->first]->Fill(*locSetIterator - *(locRFTimeSet.begin()));
340  }
341 */
342 
343  //MAKE/FILL ROC HISTOGRAMS
344  japp->RootWriteLock(); //ACQUIRE ROOT LOCK!!
345  {
346  //roc info consistency:
347  //compare roc infos: delta-t of each to avg-t-exlcuding-itself //separate histograms for each delta
348  for(size_t loc_i = 0; loc_i < locCODAROCInfos.size(); ++loc_i)
349  {
350  //double has precision of 14 digits: 4ns per tick, 14 digits = 111 hours, way longer than a run
351  double locTotalROCTimeStamp_AllButThis = double(locTotalROCTimeStamp - locCODAROCInfos[loc_i]->timestamp);
352  double locDeltaT = 4.0*(double(locCODAROCInfos[loc_i]->timestamp) - locTotalROCTimeStamp_AllButThis/(double(locCODAROCInfos.size() - 1)));
353 
354  //if histogram doesn't exist yet, create it
355  uint32_t locROCID = locCODAROCInfos[loc_i]->rocid;
356  if(dHistMap_ROCInfoDeltaT.find(locROCID) == dHistMap_ROCInfoDeltaT.end())
357  {
358  dROCTIDirectory->cd();
359  ostringstream locROCIDStream;
360  locROCIDStream << locROCID;
361  string locHistName = string("ROCDeltaTs_ROC") + locROCIDStream.str();
362  string locHistTitle = string(";#Deltat (ROC ") + locROCIDStream.str() + string(" - #mu_{Other ROCs}) (ns)");
363  dHistMap_ROCInfoDeltaT[locROCID] = new TH1I(locHistName.c_str(), locHistTitle.c_str(), 33, -16.5, 16.5);
364  gDirectory->Cd("/");
365  }
366 
367  //fill
368  dHistMap_ROCInfoDeltaT[locROCID]->Fill(locDeltaT);
369  }
370  }
371  japp->RootUnLock(); //RELEASE ROOT LOCK!!
372 
373  // FILL HISTOGRAMS
374  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
375  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
376  {
377  //num rf signals
378  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
379  {
380  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
381  locTDCIterator = locRFTimes.find(locSystem);
382  if(locTDCIterator == locRFTimes.end())
383  dHistMap_NumSignals[locSystem]->Fill(0);
384  else
385  dHistMap_NumSignals[locSystem]->Fill(locTDCIterator->second.size());
386  }
387 
388  //rf hits found, # missing, adjacent delta-t's
389  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
390  {
391  set<double>& locRFTimeSet = locTDCIterator->second;
392  DetectorSystem_t locSystem = locTDCIterator->first;
393  double locHitPeriod = dRFSignalPeriod*dRFSamplingFactor[locSystem];
394 
395  set<double>::iterator locSetIterator = locRFTimeSet.begin();
396  size_t locNumHitsMissing = 0;
397  for(; locSetIterator != locRFTimeSet.end(); ++locSetIterator)
398  {
399  double locTimeFromStart1 = *locSetIterator - *(locRFTimeSet.begin());
400  size_t locHitIndex1 = size_t(locTimeFromStart1 / locHitPeriod + 0.5);
401  dHistMap_RFHitsFound[locSystem]->Fill(locHitIndex1);
402 
403  int locExpectedHitIndex = locHitIndex1 + 1;
404  set<double>::iterator locOtherIterator = locSetIterator;
405  for(++locOtherIterator; locOtherIterator != locRFTimeSet.end(); ++locOtherIterator)
406  {
407  double locTimeFromStart2 = *locOtherIterator - *(locRFTimeSet.begin());
408  size_t locHitIndex2 = size_t(locTimeFromStart2 / locHitPeriod + 0.5);
409  locNumHitsMissing += (int(locHitIndex2) > locExpectedHitIndex) ? locHitIndex2 - locExpectedHitIndex : 0;
410  if(locHitIndex2 >= dMaxDeltaTHits[locSystem])
411  break; //no histogram for this pair
412 
413  pair<size_t, size_t> locTimePair(locHitIndex1, locHitIndex2);
414  double locDeltaT = *locOtherIterator - *locSetIterator;
415  dHistMap_AdjacentRFDeltaTs[locSystem][locTimePair]->Fill(locDeltaT);
416 
417  locExpectedHitIndex = locHitIndex2 + 1;
418  }
419  dHistMap_NumRFHitsMissing[locSystem]->Fill(locNumHitsMissing);
420  }
421  }
422 
423  //RF signal frequency
424  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
425  {
426  set<double>& locRFTimeSet = locTDCIterator->second;
427  if(locRFTimeSet.size() < 2)
428  continue;
429  DetectorSystem_t locSystem = locTDCIterator->first;
430  double locDeltaT = *(locRFTimeSet.rbegin()) - *(locRFTimeSet.begin());
431  double locRFSignalFrequency = locDeltaT/(double(dRFSamplingFactor[locSystem]*(locRFTimeSet.size() - 1)));
432  dHistMap_RFSignalPeriod[locSystem]->Fill(locRFSignalFrequency);
433  }
434 
435  //RF beam bunch period
436  for(size_t loc_i = 0; loc_i < locTAGHDeltaTs.size(); ++loc_i)
437  dHist_RFBeamBunchPeriod->Fill(locTAGHDeltaTs[loc_i]);
438 
439  //TDCs: RF / First RF
440  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
441  {
442  set<double>& locRFTimeSet = locTDCIterator->second;
443  set<double>::iterator locSetIterator = locRFTimeSet.begin();
444  for(++locSetIterator; locSetIterator != locRFTimeSet.end(); ++locSetIterator)
445  dHistMap_RFFirstTimeDeltaT[locTDCIterator->first]->Fill(*locSetIterator - *(locRFTimeSet.begin()));
446  }
447 
448  //TDCs / Self
449  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
450  {
451  set<double>& locRFTimeSet = locTDCIterator->second;
452  if(locRFTimeSet.size() < 2)
453  continue;
454  set<double>::iterator locSetIterator = locRFTimeSet.begin();
455  ++locSetIterator;
456  double locShiftedRFTime = dRFTimeFactory->Step_TimeToNearInputTime(*locSetIterator, *(locRFTimeSet.begin()));
457  dHistMap_SelfResolution[locTDCIterator->first]->Fill(locShiftedRFTime - *(locRFTimeSet.begin()));
458  }
459 
460  //RF / TAGH
461  for(size_t loc_i = 0; loc_i < locTAGHHits.size(); ++loc_i)
462  {
463  if(!locTAGHHits[loc_i]->has_TDC)
464  continue;
465 
466  //TDCs
467  for(locTDCIterator = locRFTimes.begin(); locTDCIterator != locRFTimes.end(); ++locTDCIterator)
468  {
469  set<double>& locRFTimeSet = locTDCIterator->second;
470  set<double>::iterator locSetIterator = locRFTimeSet.begin();
471  for(; locSetIterator != locRFTimeSet.end(); ++locSetIterator)
472  {
473  double locShiftedRFTime = dRFTimeFactory->Step_TimeToNearInputTime(*locSetIterator, locTAGHHits[loc_i]->time_tdc);
474  dHistMap_RFTaggerDeltaT[locTDCIterator->first]->Fill(locShiftedRFTime - locTAGHHits[loc_i]->time_tdc);
475  }
476  }
477  }
478 
479  //TDCs / TDCs
480  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
481  {
482  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
483  locTDCIterator = locRFTimes.find(locSystem);
484  if(locTDCIterator == locRFTimes.end())
485  continue;
486  set<double>& locRFTimeSet = locTDCIterator->second;
487 
488  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
489  {
490  DetectorSystem_t locSystem2 = dRFSignalSystems[loc_j];
491  map<DetectorSystem_t, set<double> >::iterator locTDCIterator2 = locRFTimes.find(locSystem2);
492  if(locTDCIterator2 == locRFTimes.end())
493  continue;
494 
495  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(locSystem, locSystem2);
496  if(dHistMap_RFRFDeltaTs.find(locSystemPair) == dHistMap_RFRFDeltaTs.end())
497  continue;
498  set<double>& locRFTimeSet2 = locTDCIterator2->second;
499 
500  double locShiftedRFTime = dRFTimeFactory->Step_TimeToNearInputTime(*(locRFTimeSet.begin()), *(locRFTimeSet2.begin()));
501  double locDeltaT = locShiftedRFTime - *(locRFTimeSet2.begin());
502  dHistMap_RFRFDeltaTs[locSystemPair]->Fill(locDeltaT);
503  }
504  }
505 
506  //Average TDCs / TDCs
507  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
508  {
509  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
510  locTDCIterator = locRFTimes.find(locSystem);
511  if(locTDCIterator == locRFTimes.end())
512  continue;
513 
514  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
515  {
516  DetectorSystem_t locSystem2 = dRFSignalSystems[loc_j];
517  map<DetectorSystem_t, set<double> >::iterator locTDCIterator2 = locRFTimes.find(locSystem2);
518  if(locTDCIterator2 == locRFTimes.end())
519  continue;
520 
521  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(locSystem, locSystem2);
522  if(dHistMap_AverageRFRFDeltaTs.find(locSystemPair) == dHistMap_AverageRFRFDeltaTs.end())
523  continue;
524 
525  double locShiftedAverageRFTime = dRFTimeFactory->Step_TimeToNearInputTime(locAverageRFTimes[locTDCIterator->first], locAverageRFTimes[locTDCIterator2->first]);
526  double locDeltaT = locShiftedAverageRFTime - locAverageRFTimes[locTDCIterator2->first];
527  dHistMap_AverageRFRFDeltaTs[locSystemPair]->Fill(locDeltaT);
528  }
529  }
530 
531  //Absolute TDCs / TDCs
532  for(size_t loc_i = 0; loc_i < dRFSignalSystems.size(); ++loc_i)
533  {
534  DetectorSystem_t locSystem = dRFSignalSystems[loc_i];
535  locTDCIterator = locRFTimes.find(locSystem);
536  if(locTDCIterator == locRFTimes.end())
537  continue;
538  set<double>& locRFTimeSet = locTDCIterator->second;
539 
540  for(size_t loc_j = loc_i + 1; loc_j < dRFSignalSystems.size(); ++loc_j)
541  {
542  DetectorSystem_t locSystem2 = dRFSignalSystems[loc_j];
543  map<DetectorSystem_t, set<double> >::iterator locTDCIterator2 = locRFTimes.find(locSystem2);
544  if(locTDCIterator2 == locRFTimes.end())
545  continue;
546 
547  pair<DetectorSystem_t, DetectorSystem_t> locSystemPair(locSystem, locSystem2);
548  if(dHistMap_AbsoluteRFRFDeltaTs.find(locSystemPair) == dHistMap_AbsoluteRFRFDeltaTs.end())
549  continue;
550  set<double>& locRFTimeSet2 = locTDCIterator2->second;
551 
552  double locDeltaT = *(locRFTimeSet.begin()) - *(locRFTimeSet2.begin());
553  dHistMap_AbsoluteRFRFDeltaTs[locSystemPair]->Fill(locDeltaT);
554  }
555  }
556  }
557  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
558 
559  return NOERROR;
560 }
561 
563 {
564  // This is called whenever the run number changes, before it is
565  // changed to give you a chance to clean up before processing
566  // events from the next run number.
567  return NOERROR;
568 }
569 
571 {
572  // Called before program exit after event processing is finished.
573  return NOERROR;
574 }
jerror_t evnt(jana::JEventLoop *locEventLoop, uint64_t eventnumber)
Called every event.
map< pair< DetectorSystem_t, DetectorSystem_t >, TH1I * > dHistMap_RFRFDeltaTs
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
char string[256]
map< DetectorSystem_t, double > dRFSamplingFactor
axes Fill(100, 100)
map< DetectorSystem_t, TH1I * > dHistMap_RFHitsFound
Definition: GlueX.h:29
double Convert_TDCToTime(const DRFTDCDigiTime *locRFTDCDigiTime, const DTTabUtilities *locTTabUtilities) const
map< DetectorSystem_t, TH1I * > dHistMap_RFFirstTimeDeltaT
vector< DetectorSystem_t > dRFSignalSystems
map< DetectorSystem_t, TH1I * > dHistMap_SelfResolution
DetectorSystem_t
Definition: GlueX.h:15
map< DetectorSystem_t, TH1I * > dHistMap_RFSignalPeriod
map< DetectorSystem_t, TH1I * > dHistMap_RFTaggerDeltaT
JApplication * japp
map< DetectorSystem_t, size_t > dMaxDeltaTHits
map< uint32_t, TH1I * > dHistMap_ROCInfoDeltaT
InitPlugin_t InitPlugin
Definition: GlueX.h:20
Definition: GlueX.h:18
map< DetectorSystem_t, TH1I * > dHistMap_NumSignals
jerror_t brun(jana::JEventLoop *locEventLoop, int32_t runnumber)
Called everytime a new run number is detected.
map< DetectorSystem_t, map< pair< size_t, size_t >, TH1I * > > dHistMap_AdjacentRFDeltaTs
const char * SystemName(DetectorSystem_t sys)
Definition: GlueX.h:38
jerror_t init(void)
Called once at program start.
map< pair< DetectorSystem_t, DetectorSystem_t >, TH1I * > dHistMap_AbsoluteRFRFDeltaTs
double Step_TimeToNearInputTime(double locTimeToStep, double locTimeToStepTo) const
DRFTime_factory * dRFTimeFactory
map< pair< DetectorSystem_t, DetectorSystem_t >, TH1I * > dHistMap_AverageRFRFDeltaTs
Definition: GlueX.h:32
map< DetectorSystem_t, TH1I * > dHistMap_NumRFHitsMissing
jerror_t fini(void)
Called after last event of last event source has been processed.