Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_TAGH_online.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_FCAL_online.cc
4 // Created: Fri Nov 9 11:58:09 EST 2012
5 // Creator: wolin (on Linux stan.jlab.org 2.6.32-279.11.1.el6.x86_64 x86_64)
6 
7 
8 #include <stdint.h>
9 #include <vector>
10 
12 #include <JANA/JApplication.h>
13 
14 using namespace std;
15 using namespace jana;
16 
17 #include "TTAB/DTTabUtilities.h"
18 #include <TAGGER/DTAGHTDCDigiHit.h>
19 #include <TAGGER/DTAGHDigiHit.h>
20 #include <TAGGER/DTAGHHit.h>
21 #include <TAGGER/DTAGHGeometry.h>
22 #include <DAQ/Df250PulseData.h>
23 #include <DAQ/Df250WindowRawData.h>
24 #include <DAQ/DEPICSvalue.h>
25 
26 #include <TDirectory.h>
27 #include <TH1.h>
28 #include <TH2.h>
29 #include <TProfile.h>
30 
31 // Nslots: total number of TAGH counter slots
32 // 274 total (218 have counters for default GlueX configuration)
34 // counter (slot) id = HV id for hv id from 1-131 (upstream of microscope)
35 // counter (slot) id = 41 + HV id for hv id from 132-233 (downstream of microscope)
36 const int NHVchannels = 233;
37 const int NupstreamCounters = 131;
38 //
39 const int NmultBins = 300; // number of bins for multiplicity histograms
40 //
41 double beam_current = -1.0;
42 const double counts_cut = 1000.0;
43 // root hist pointers
44 static TH1I *tagh_num_events;
45 static TH1I *hBeamCurrent;
46 static TH2F *hHit_HasTDCvsHasADC;
47 static TH1I *hHit_RawNHits;
48 static TH1I *hHit_NHits;
49 static TH1I *hHit_NHits_us;
50 static TH1I *hHit_NHits_ds;
51 static TH1F *hHit_Occupancy;
52 static TH2I *hHit_HVidVsSlotID;
53 static TH1F *hHit_Energy;
54 static TH2I *hHit_EnergyVsSlotID;
55 static TH1I *hHit_Integral;
56 static TH2I *hHit_IntegralVsSlotID;
57 //static TH1I *hHit_fadcNpe;
58 static TH1I *hHit_fadcTime;
59 static TH2I *hHit_fadcTimeVsSlotID;
60 static TH1I *hHit_Time;
61 static TH2F *hHit_TimeVsSlotID;
62 static TH2F *hHit_TimeVsEnergy;
63 static TH2I *hHit_TimeVsIntegral;
64 static TH1I *hHit_tdcTime;
65 static TH2I *hHit_tdcTimeVsSlotID;
68 //
69 static TH1I *hDigiHit_NfadcHits;
73 static TH1I *hDigiHit_Pedestal;
74 static TProfile *hDigiHit_PedestalVsSlotID;
76 static TH1I *hDigiHit_PulseNumber;
79 static TH1I *hDigiHit_RawPeak;
81 static TH1I *hDigiHit_RawIntegral;
84 static TH2I *hDigiHit_PeakVsSlotID;
87 static TH1I *hDigiHit_PulseTime;
88 static TH1I *hDigiHit_fadcTime;
91 static TH1I *hDigiHit_NtdcHits;
93 static TH1I *hDigiHit_tdcOccupancy;
94 static TH1I *hDigiHit_tdcRawTime;
95 static TH1I *hDigiHit_tdcTime;
101 //
108 //----------------------------------------------------------------------------------
109 
110 
111 // Routine used to create our JEventProcessor
112 extern "C"{
113  void InitPlugin(JApplication *app){
114  InitJANAPlugin(app);
115  app->AddProcessor(new JEventProcessor_TAGH_online());
116  }
117 }
118 
119 
120 //----------------------------------------------------------------------------------
121 
122 
124 }
125 
126 
127 //----------------------------------------------------------------------------------
128 
129 
131 }
132 
133 
134 //----------------------------------------------------------------------------------
135 
137 
138  // create root folder for tagh and cd to it, store main dir
139  TDirectory *mainDir = gDirectory;
140  TDirectory *taghDir = gDirectory->mkdir("TAGH");
141  taghDir->cd();
142 
143  // book hists
144  tagh_num_events = new TH1I("tagh_num_events","TAGH number of events",1,0.5,1.5);
145  hBeamCurrent = new TH1I("BeamCurrent","Beam current;beam current [nA]",600,0.0,300.0);
146  // hit-level hists (after calibration)
147  gDirectory->mkdir("Hit")->cd();
148  hHit_HasTDCvsHasADC = new TH2F("Hit_HasTDCvsHasADC","TAGH has TDC? vs. has ADC?;fADC status;TDC status",2,-0.5,1.5,2,-0.5,1.5);
149  hHit_RawNHits = new TH1I("Hit_RawNHits","TAGH hit multiplicity (raw);hits (raw);events",NmultBins,0.5,0.5+NmultBins);
150  //
151  hHit_NHits = new TH1I("Hit_NHits","TAGH hit multiplicity;hits;events",NmultBins,0.5,0.5+NmultBins);
152  hHit_NHits_us = new TH1I("Hit_NHits_us","TAGH upstream hit multiplicity;upstream hits;events",NmultBins,0.5,0.5+NmultBins);
153  hHit_NHits_ds = new TH1I("Hit_NHits_ds","TAGH downstream hit multiplicity;downstream hits;events",NmultBins,0.5,0.5+NmultBins);
154  hHit_Occupancy = new TH1F("Hit_Occupancy","TAGH occupancy;counter (slot) ID;hits / counter",Nslots,0.5,0.5+Nslots);
155  hHit_HVidVsSlotID = new TH2I("Hit_HVidVsSlotID","TAGH high voltage channel ID vs. counter ID;counter (slot) ID;high voltage channel ID",Nslots,0.5,0.5+Nslots,NHVchannels,0.5,0.5+NHVchannels);
156  hHit_Energy = new TH1F("Hit_Energy","TAGH energy;photon energy [GeV];hits / counter",120,0.0,12.0);
157  hHit_EnergyVsSlotID = new TH2I("Hit_EnergyVsSlotID","TAGH energy vs. counter ID;counter (slot) ID;photon energy [GeV]",Nslots,0.5,0.5+Nslots,120,0.0,12.0);
158  hHit_Integral = new TH1I("Hit_Integral","TAGH fADC pulse integral;pulse integral;hits",1000,0.0,30000.0);
159  hHit_IntegralVsSlotID = new TH2I("Hit_IntegralVsSlotID","TAGH fADC pulse integral vs. counter ID;counter (slot) ID;pulse integral",Nslots,0.5,0.5+Nslots,1000,0.0,30000.0);
160  //hHit_fadcNpe = new TH1I("Hit_fadcNpe","TAGH fADC number of photoelectrons;photoelectrons;hits",1000,0.0,30000.0);
161  hHit_fadcTime = new TH1I("Hit_fadcTime","TAGH fADC time;time [ns];hits / 400 ps",2000,-400.0,400.0);
162  hHit_fadcTimeVsSlotID = new TH2I("Hit_fadcTimeVsSlotID","TAGH fADC time vs. counter ID;counter (slot) ID;time [ns]",Nslots,0.5,0.5+Nslots,2000,-400.0,400.0);
163  hHit_Time = new TH1I("Hit_Time","TAGH time;time [ns];hits / 400 ps",2000,-400.0,400.0);
164  hHit_TimeVsSlotID = new TH2F("Hit_TimeVsSlotID","TAGH time vs. counter ID;counter (slot) ID;time [ns]",Nslots,0.5,0.5+Nslots,2000,-400.0,400.0);
165  hHit_TimeVsEnergy = new TH2F("Hit_TimeVsEnergy","TAGH time vs. energy;energy [GeV];time [ns]",120,0.0,12.0,2000,-400.0,400.0);
166  hHit_TimeVsIntegral = new TH2I("Hit_TimeVsIntegral","TAGH time vs. integral;pulse integral;time [ns]",500,0.0,15000.0,2000,-400.0,400.0);
167  hHit_tdcTime = new TH1I("Hit_tdcTime","TAGH TDC time;time [ns];hits / 400 ps",2000,-400.0,400.0);
168  hHit_tdcTimeVsSlotID = new TH2I("Hit_tdcTimeVsSlotID","TAGH TDC time vs. counter ID;counter (slot) ID;time [ns]",Nslots,0.5,0.5+Nslots,2000,-400.0,400.0);
169  hHit_tdcadcTimeDiffVsSlotID = new TH2I("Hit_tdcadcTimeDiffVsSlotID","TAGH TDC/ADC time difference vs. counter ID;counter (slot) ID;time(TDC) - time(ADC) [ns]",Nslots,0.5,0.5+Nslots,200,-40.0,40.0);
170  hHit_tdcadcTimeDiffVsIntegral = new TH2I("Hit_tdcadcTimeDiffVsIntegral","TAGH TDC/ADC time difference vs. integral;pulse integral;time(TDC) - time(ADC) [ns]",500,0.0,15000.0,200,-40.0,40.0);
171  // digihit-level hists
172  taghDir->cd();
173  gDirectory->mkdir("DigiHit")->cd();
174  hDigiHit_NfadcHits = new TH1I("DigiHit_NfadcHits","TAGH fADC hit multiplicity;raw hits;events",NmultBins,0.5,0.5+NmultBins);
175  hDigiHit_NSamplesPedestal = new TH1I("DigiHit_NSamplesPedestal","TAGH fADC pedestal samples;pedestal samples;raw hits",50,-0.5,49.5);
176  hDigiHit_Pedestal = new TH1I("DigiHit_Pedestal","TAGH fADC pedestals;pedestal [fADC counts];raw hits",200,0.0,200.0);
177  hDigiHit_PedestalVsSlotID = new TProfile("DigiHit_PedestalVsSlotID","TAGH pedestal vs. counter ID;counter (slot) ID;average pedestal [fADC counts]",Nslots,0.5,0.5+Nslots,"s");
178  hDigiHit_QualityFactor = new TH1I("DigiHit_QualityFactor","TAGH fADC quality factor;quality factor;raw hits",4,-0.5,3.5);
179  hDigiHit_PulseNumber = new TH1I("DigiHit_PulseNumber","TAGH fADC pulse number;pulse number;raw hits",4,-0.5,3.5);
180  hDigiHit_PulseNumberVsSlotID = new TH2I("DigiHit_PulseNumberVsSlotID","TAGH fADC pulse number vs. counter ID;counter (slot) ID;pulse number;raw hits",Nslots,0.5,0.5+Nslots,4,-0.5,3.5);
181  hDigiHit_fadcOccupancy = new TH1I("DigiHit_fadcOccupancy","TAGH fADC hit occupancy;counter (slot) ID;raw hits / counter",Nslots,0.5,0.5+Nslots);
182  hDigiHit_RawPeak = new TH1I("DigiHit_RawPeak","TAGH fADC pulse peak (raw);pulse peak (raw);raw hits",410,0.0,4100.0);
183  hDigiHit_RawPeakVsSlotID = new TH2I("DigiHit_RawPeakVsSlotID","TAGH fADC pulse peak (raw) vs. counter ID;counter (slot) ID;pulse peak (raw)",Nslots,0.5,0.5+Nslots,410,0.0,4100.0);
184  hDigiHit_RawIntegral = new TH1I("DigiHit_RawIntegral","TAGH fADC pulse integral (raw);pulse integral (raw);raw hits",1000,0.0,30000.0);
185  hDigiHit_RawIntegralVsSlotID = new TH2I("DigiHit_RawIntegralVsSlotID","TAGH fADC pulse integral (raw) vs. counter ID;counter (slot) ID;pulse integral (raw)",Nslots,0.5,0.5+Nslots,1000,0.0,30000.0);
186  hDigiHit_NSamplesIntegral = new TH1I("DigiHit_NSamplesIntegral","TAGH fADC integral samples;integral samples;raw hits",60,-0.5,59.5);
187  hDigiHit_PeakVsSlotID = new TH2I("DigiHit_PeakVsSlotID","TAGH fADC pulse peak vs. counter ID;counter (slot) ID;pulse peak",Nslots,0.5,0.5+Nslots,410,0.0,4100.0);
188  hDigiHit_IntegralVsPeak = new TH2I("DigiHit_IntegralVsPeak","TAGH fADC pulse integral vs. peak;pulse peak;pulse integral",410,0.0,4100.0,1000,0.0,30000.0);
189  hDigiHit_IntegralVsSlotID = new TH2I("DigiHit_IntegralVsSlotID","TAGH fADC pulse integral vs. counter ID;counter (slot) ID;pulse integral",Nslots,0.5,0.5+Nslots,1000,0.0,30000.0);
190  hDigiHit_PulseTime = new TH1I("DigiHit_PulseTime","TAGH fADC pulse time;pulse time [62.5 ps];raw hits",1000,0.0,6500.0);
191  hDigiHit_fadcTime = new TH1I("DigiHit_fadcTime","TAGH fADC pulse time;pulse time [ns];raw hits / 400 ps",1000,0.0,400.0);
192  hDigiHit_fadcTimeVsSlotID = new TH2I("DigiHit_fadcTimeVsSlotID","TAGH fADC pulse time vs. counter ID;counter (slot) ID;pulse time [ns]",Nslots,0.5,0.5+Nslots,1000,0.0,400.0);
193  hDigiHit_fadcTimeVsIntegral = new TH2I("DigiHit_fadcTimeVsIntegral","TAGH fADC pulse time vs. integral;pulse integral;pulse time [ns]",500,0.0,15000.0,1000,0.0,400.0);
194  hDigiHit_NtdcHits = new TH1I("DigiHit_NtdcHits","TAGH TDC hit multiplicity;raw hits;events",NmultBins,0.5,0.5+NmultBins);
195  hDigiHit_NtdcHitsVsSlotID = new TH2I("DigiHit_NtdcHitsVsSlotID","TAGH TDC hit multiplicity vs. counter ID;counter ID;raw hits",Nslots,0.5,0.5+Nslots,8,0.5,8.5);
196  hDigiHit_tdcOccupancy = new TH1I("DigiHit_tdcOccupancy","TAGH TDC hit occupancy;counter (slot) ID;raw hits / counter",Nslots,0.5,0.5+Nslots);
197  hDigiHit_tdcRawTime = new TH1I("DigiHit_tdcRawTime","TAGH TDC raw time;time [60 ps];raw hits",1000,0.0,65500.0);
198  hDigiHit_tdcTime = new TH1I("DigiHit_tdcTime","TAGH TDC time;time [ns];raw hits / 400 ps",2000,0.0,800.0);
199  hDigiHit_tdcTimeVsSlotID = new TH2I("DigiHit_tdcTimeVsSlotID","TAGH TDC time vs. counter ID;counter ID;TDC time [ns]",Nslots,0.5,0.5+Nslots,2000,0.0,800.0);
200  hDigiHit_tdcTimeVsfadcTime = new TH2I("DigiHit_tdcTimeVsfadcTime","TAGH TDC time vs. ADC time;fADC time [ns];TDC time [ns]",400,0.0,400.0,800,0.0,800.0);
201  hDigiHit_tdcadcTimeDiff = new TH1I("DigiHit_tdcadcTimeDiff","TAGH TDC/ADC time difference;time(TDC) - time(ADC) [ns];raw hits / 400 ps",1000,0.0,400.0);
202  hDigiHit_tdcadcTimeDiffVsSlotID = new TH2I("DigiHit_tdcadcTimeDiffVsSlotID","TAGH TDC/ADC time difference vs. counter ID;counter ID;time(TDC) - time(ADC) [ns]",Nslots,0.5,0.5+Nslots,1000,0.0,400.0);
203  hDigiHit_tdcadcTimeDiffVsIntegral = new TH2I("DigiHit_tdcadcTimeDiffVsIntegral","TAGH TDC/ADC time difference vs. pulse integral;pulse integral;time(TDC) - time(ADC) [ns]",500,0.0,15000.0,1000,0.0,400.0);
204  // allow multiple peaks in window and cut on ADC counts
205  hDigiHit_NfadcHits_multiPeak = new TH1I("DigiHit_NfadcHits_multiPeak","TAGH fADC hit multiplicity (> 400 counts,multiple peaks allowed);raw hits;events",NmultBins,0.5,0.5+NmultBins);
206  hDigiHit_NfadcHitsVsSlotID = new TH2I("DigiHit_NfadcHitsVsSlotID","TAGH fADC hit multiplicity vs. counter ID (> 400 counts);counter ID;raw hits",Nslots,0.5,0.5+Nslots,8,0.5,8.5);
207  // digihit-level hists after cut on ADC integral counts
208  hDigiHit_NfadcHits_cut = new TH1I("DigiHit_NfadcHits_cut","TAGH fADC hit multiplicity (> 1k ADC integral counts);raw hits;events",NmultBins,0.5,0.5+NmultBins);
209  hDigiHit_fadcOccupancy_cut = new TH1I("DigiHit_fadcOccupancy_cut","TAGH fADC hit occupancy (> 1k ADC integral counts);counter (slot) ID;raw hits / counter",Nslots,0.5,0.5+Nslots);
210  hDigiHit_fadcTime_cut = new TH1I("DigiHit_fadcTime_cut","TAGH fADC pulse time (> 1k ADC integral counts);pulse time [ns];raw hits / 400 ps",1000,0.0,400.0);
211  hDigiHit_fadcTimeVsSlotID_cut = new TH2I("DigiHit_fadcTimeVsSlotID_cut","TAGH fADC pulse time vs. counter ID (> 1k ADC integral counts);counter (slot) ID;pulse time [ns]",Nslots,0.5,0.5+Nslots,1000,0.0,400.0);
212  hDigiHit_fadcTimeVsQF_cut = new TH2I("DigiHit_fadcTimeVsQF_cut","TAGH fADC pulse time vs. quality factor (> 1k ADC integral counts);fADC quality factor;pulse time [ns]",4,-0.5,3.5,1000,0.0,400.0);
213  hDigiHit_adctdcMatchesVsSlotID_cut = new TH2I("DigiHit_adctdcMatchesVsSlotID_cut","TAGH #TDC matches / fADC hit vs. counter ID (> 1k ADC integral counts);counter ID;#TDC matches / fADC hit",Nslots,0.5,0.5+Nslots,8,-0.5,7.5);
214  // back to main dir
215  mainDir->cd();
216 
217  return NOERROR;
218 }
219 
220 //----------------------------------------------------------------------------------
221 
222 
223 jerror_t JEventProcessor_TAGH_online::brun(JEventLoop *eventLoop, int32_t runnumber) {
224  // This is called whenever the run number changes
225  // extract the TAGH geometry
226  vector<const DTAGHGeometry*> taghGeomVect;
227  eventLoop->Get(taghGeomVect);
228  if (taghGeomVect.size() == 0) return OBJECT_NOT_AVAILABLE;
229  const DTAGHGeometry& taghGeom = *(taghGeomVect[0]);
230  // get photon energy bin low of each counter for energy histogram binning
231  double Elows[Nslots + 1];
232  double slots[Nslots + 1];
233  for (int i = 0; i < Nslots; i++) {
234  Elows[i] = taghGeom.getElow(Nslots - i);
235  slots[i] = 0.5 + i;
236  }
237  // add the upper limit
238  Elows[Nslots] = taghGeom.getEhigh(1);
239  slots[Nslots] = 0.5 + Nslots;
240 
241  const int Ntime = 2000;
242  double Tlows[Ntime + 1];
243  for (int i = 0; i <= Ntime; i++) {
244  Tlows[i] = -400.0 + i*0.4;
245  }
246 
247  // FILL HISTOGRAMS
248  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
249  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
250 
251  if (hHit_Energy->GetEntries() == 0) hHit_Energy->SetBins(Nslots,Elows);
252  if (hHit_EnergyVsSlotID->GetEntries() == 0) hHit_EnergyVsSlotID->SetBins(Nslots,slots,Nslots,Elows);
253  if (hHit_TimeVsEnergy->GetEntries() == 0) hHit_TimeVsEnergy->SetBins(Nslots,Elows,Ntime,Tlows);
254 
255  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
256 
257  return NOERROR;
258 }
259 
260 
261 //----------------------------------------------------------------------------------
262 
263 
264 jerror_t JEventProcessor_TAGH_online::evnt(JEventLoop *eventLoop, uint64_t eventnumber) {
265  // This is called for every event. Use of common resources like writing
266  // to a file or filling a histogram should be mutex protected. Using
267  // loop-Get(...) to get reconstructed objects (and thereby activating the
268  // reconstruction algorithm) should be done outside of any mutex lock
269  // since multiple threads may call this method at the same time.
270 
271  // Get TAGH hits and digihits
272  vector<const DTAGHHit*> hits;
273  eventLoop->Get(hits, "Calib");
274  vector<const DTAGHDigiHit*> digihits;
275  eventLoop->Get(digihits);
276  vector<const DTAGHTDCDigiHit*> tdcdigihits;
277  eventLoop->Get(tdcdigihits);
278 
279  const DTTabUtilities* ttabUtilities = nullptr;
280  eventLoop->GetSingle(ttabUtilities);
281 
282  // Cache pulse data and window raw data objects
283  map< const DTAGHDigiHit*, pair<const Df250PulseData*, const Df250WindowRawData*> > pd_wrd_cache;
284  for (const auto& hit : digihits) {
285  const Df250PulseData* pd = nullptr;
286  const Df250WindowRawData* wrd = nullptr;
287  hit->GetSingle(pd);
288  if (pd != nullptr) pd->GetSingle(wrd);
289  pd_wrd_cache[hit] = pair<const Df250PulseData*, const Df250WindowRawData*>(pd, wrd);
290  }
291 
292  // Get beam current
293  vector<const DEPICSvalue*> epicsVals;
294  eventLoop->Get(epicsVals);
295  for (const auto& ev : epicsVals) {
296  if (ev->name == "IBCAD00CRCUR6") {
297  beam_current = ev->fval;
298  break;
299  }
300  }
301 
302  // FILL HISTOGRAMS
303  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
304  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
305 
306  hBeamCurrent->Fill(beam_current);
307  if (digihits.size() > 0 || tdcdigihits.size() > 0)
308  tagh_num_events->Fill(1);
309 
310  hHit_RawNHits->Fill(hits.size());
311  hDigiHit_NfadcHits->Fill(digihits.size());
312  hDigiHit_NtdcHits->Fill(tdcdigihits.size());
313  int NfadcHits_cut = 0;
314  int NHits_hasADC = 0; int NHits_hasADC_us = 0; int NHits_hasADC_ds = 0;
315  int Nadc[Nslots];
316  for (int i = 0; i < Nslots; i++) Nadc[i] = 0;
317  for (const auto& hit : digihits) {
318  double ped = (double)hit->pedestal/hit->nsamples_pedestal;
319  hDigiHit_NSamplesPedestal->Fill(hit->nsamples_pedestal);
320  hDigiHit_Pedestal->Fill(ped);
321  hDigiHit_RawPeak->Fill(hit->pulse_peak);
322  if (ped == 0.0 || hit->pulse_peak == 0) continue;
323  hDigiHit_PedestalVsSlotID->Fill(hit->counter_id,ped);
324  hDigiHit_fadcOccupancy->Fill(hit->counter_id);
325  hDigiHit_RawPeakVsSlotID->Fill(hit->counter_id,hit->pulse_peak);
326  hDigiHit_RawIntegral->Fill(hit->pulse_integral);
327  hDigiHit_RawIntegralVsSlotID->Fill(hit->counter_id,hit->pulse_integral);
328  hDigiHit_NSamplesIntegral->Fill(hit->nsamples_integral);
329  double pI = hit->pulse_integral-hit->nsamples_integral*ped;
330  hDigiHit_IntegralVsSlotID->Fill(hit->counter_id,pI);
331  hDigiHit_IntegralVsPeak->Fill(hit->pulse_peak-ped,pI);
332  hDigiHit_PeakVsSlotID->Fill(hit->counter_id,hit->pulse_peak-ped);
333  hDigiHit_PulseTime->Fill(hit->pulse_time);
334  double t_ns = 0.0625*hit->pulse_time;
335  hDigiHit_fadcTime->Fill(t_ns);
336  hDigiHit_fadcTimeVsSlotID->Fill(hit->counter_id,t_ns);
337  hDigiHit_fadcTimeVsIntegral->Fill(pI,t_ns);
338  hDigiHit_QualityFactor->Fill(hit->QF);
339  const Df250PulseData* pd = pd_wrd_cache[hit].first;
340  int pN = (pd != nullptr) ? pd->pulse_number : -1;
341  if (pN > -1) hDigiHit_PulseNumber->Fill(pN);
342  if (pN > -1) hDigiHit_PulseNumberVsSlotID->Fill(hit->counter_id,pN);
343  if (pI > counts_cut) {
344  NfadcHits_cut++;
345  hDigiHit_fadcOccupancy_cut->Fill(hit->counter_id);
346  hDigiHit_fadcTime_cut->Fill(t_ns);
347  hDigiHit_fadcTimeVsSlotID_cut->Fill(hit->counter_id,t_ns);
348  hDigiHit_fadcTimeVsQF_cut->Fill(hit->QF,t_ns);
349  }
350  const Df250WindowRawData* wrd = pd_wrd_cache[hit].second;
351  size_t threshold = 400;
352  if (wrd != nullptr) {
353  size_t prevSample = 100;
354  int Ncrosses = 0;
355  if (!wrd->invalid_samples && !wrd->overflow) {
356  for (size_t sample : wrd->samples) {
357  if (sample > threshold && prevSample < threshold) {
358  Ncrosses++;
359  }
360  prevSample = sample;
361  }
362  }
363  Nadc[hit->counter_id-1] = Ncrosses;
364  }
365  else {
366  if (hit->pulse_peak > threshold) Nadc[hit->counter_id-1]++;
367  }
368  }
369 
370  int NmultiPeak = 0;
371  if (digihits.size() > 0) {
372  for (int i = 0; i < Nslots; i++) {
373  NmultiPeak += Nadc[i];
374  hDigiHit_NfadcHitsVsSlotID->Fill(i+1,Nadc[i]);
375  }
376  }
377  hDigiHit_NfadcHits_multiPeak->Fill(NmultiPeak);
378  hDigiHit_NfadcHits_cut->Fill(NfadcHits_cut);
379  for (const auto& hit : digihits) {
380  double ped = (double)hit->pedestal/hit->nsamples_pedestal;
381  double pI = hit->pulse_integral-ped*hit->nsamples_integral;
382  if (hit->pedestal > 0 && pI > counts_cut) {
383  int matches = 0;
384  for (const auto& tdchit : tdcdigihits) {
385  if (hit->counter_id == tdchit->counter_id) {
386  matches++;
387  double T_tdc = ttabUtilities->Convert_DigiTimeToNs_F1TDC(tdchit);
388  double T_adc = 0.0625*hit->pulse_time;
389  hDigiHit_tdcTimeVsfadcTime->Fill(T_adc,T_tdc);
390  hDigiHit_tdcadcTimeDiff->Fill(T_tdc-T_adc);
391  hDigiHit_tdcadcTimeDiffVsSlotID->Fill(hit->counter_id,T_tdc-T_adc);
392  hDigiHit_tdcadcTimeDiffVsIntegral->Fill(pI,T_tdc-T_adc);
393  }
394  }
395  hDigiHit_adctdcMatchesVsSlotID_cut->Fill(hit->counter_id,matches);
396  }
397  }
398 
399  int Ntdc[Nslots];
400  for (int i = 0; i < Nslots; i++) Ntdc[i] = 0;
401  for (const auto& hit : tdcdigihits) {
402  Ntdc[hit->counter_id-1]++;
403  hDigiHit_tdcOccupancy->Fill(hit->counter_id);
404  hDigiHit_tdcRawTime->Fill(hit->time);
405  double T_tdc = ttabUtilities->Convert_DigiTimeToNs_F1TDC(hit);
406  hDigiHit_tdcTime->Fill(T_tdc);
407  hDigiHit_tdcTimeVsSlotID->Fill(hit->counter_id,T_tdc);
408  }
409  if (tdcdigihits.size() > 0) {
410  for (int i = 0; i < Nslots; i++) {
411  hDigiHit_NtdcHitsVsSlotID->Fill(i+1,Ntdc[i]);
412  }
413  }
414 
415  // Fill histograms with calibrated-hit data
416  for (const auto& hit : hits) {
417  hHit_HasTDCvsHasADC->Fill(hit->has_fADC,hit->has_TDC);
418  if (hit->has_fADC) {
419  NHits_hasADC++;
420  if (hit->counter_id<=NupstreamCounters) NHits_hasADC_us++;
421  else NHits_hasADC_ds++;
422  hHit_Occupancy->Fill(hit->counter_id);
423  int HVid = (hit->counter_id <= NupstreamCounters) ? hit->counter_id : (hit->counter_id - Nslots + NHVchannels);
424  hHit_HVidVsSlotID->Fill(hit->counter_id,HVid);
425  hHit_Energy->Fill(hit->E);
426  hHit_EnergyVsSlotID->Fill(hit->counter_id,hit->E);
427  hHit_Integral->Fill(hit->integral);
428  hHit_IntegralVsSlotID->Fill(hit->counter_id,hit->integral);
429  //hHit_fadcNpe->Fill(hit->npe_fadc);
430  hHit_fadcTime->Fill(hit->time_fadc);
431  hHit_fadcTimeVsSlotID->Fill(hit->counter_id,hit->time_fadc);
432  hHit_Time->Fill(hit->t);
433  hHit_TimeVsSlotID->Fill(hit->counter_id,hit->t);
434  hHit_TimeVsEnergy->Fill(hit->E,hit->t);
435  hHit_TimeVsIntegral->Fill(hit->integral,hit->t);
436  if (hit->has_TDC) {
437  hHit_tdcTime->Fill(hit->time_tdc);
438  hHit_tdcTimeVsSlotID->Fill(hit->counter_id,hit->time_tdc);
439  hHit_tdcadcTimeDiffVsSlotID->Fill(hit->counter_id,hit->time_tdc-hit->time_fadc);
440  hHit_tdcadcTimeDiffVsIntegral->Fill(hit->integral,hit->time_tdc-hit->time_fadc);
441  }
442  }
443  }
444  hHit_NHits->Fill(NHits_hasADC);
445  hHit_NHits_us->Fill(NHits_hasADC_us);
446  hHit_NHits_ds->Fill(NHits_hasADC_ds);
447 
448  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
449 
450  return NOERROR;
451 }
452 
453 
454 //----------------------------------------------------------------------------------
455 
456 
458  // This is called whenever the run number changes, before it is
459  // changed to give you a chance to clean up before processing
460  // events from the next run number.
461  return NOERROR;
462 }
463 
464 
465 //----------------------------------------------------------------------------------
466 
467 
469  // Called before program exit after event processing is finished.
470  return NOERROR;
471 }
472 
473 
474 //----------------------------------------------------------------------------------
475 //----------------------------------------------------------------------------------
static TH1I * hHit_NHits_us
double counts_cut
double Convert_DigiTimeToNs_F1TDC(const JObject *locTDCDigiHit) const
double getElow(unsigned int counter) const
static TH1I * hDigiHit_NSamplesIntegral
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
static TH1I * hDigiHit_NSamplesPedestal
static TH1I * hDigiHit_NfadcHits
if(locHist_BCALShowerPhiVsZ!=NULL)
static TH2I * hDigiHit_adctdcMatchesVsSlotID_cut
static TH2I * hDigiHit_IntegralVsSlotID
static TH1I * hHit_Time
static TH1I * hHit_Integral
static TH2I * hDigiHit_RawIntegralVsSlotID
static TH2I * hDigiHit_tdcadcTimeDiffVsSlotID
static TH2F * hHit_TimeVsEnergy
static TH2I * hDigiHit_NtdcHitsVsSlotID
static TH1I * tagh_num_events
const int NmultBins
vector< uint16_t > samples
const int NHVchannels
TDirectory * mainDir
Definition: p2k_hists.C:2
static TH1I * hDigiHit_fadcOccupancy
static TH2I * hHit_IntegralVsSlotID
static TH2I * hDigiHit_fadcTimeVsSlotID_cut
static TH1I * hDigiHit_NtdcHits
static TH1I * hBeamCurrent
static TH2I * hDigiHit_fadcTimeVsIntegral
JApplication * japp
static TH2I * hDigiHit_tdcadcTimeDiffVsIntegral
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
static TH1F * hHit_Occupancy
static TH1I * hDigiHit_fadcTime_cut
static TH1I * hDigiHit_fadcOccupancy_cut
static TH1I * hHit_tdcTime
static TH2I * hHit_EnergyVsSlotID
static TH1I * hDigiHit_tdcOccupancy
static TH2I * hDigiHit_IntegralVsPeak
static TH2I * hDigiHit_PeakVsSlotID
InitPlugin_t InitPlugin
static TH1F * hHit_Energy
static TH1I * hHit_fadcTime
static TH2I * hHit_TimeVsIntegral
static TH2I * hHit_tdcadcTimeDiffVsSlotID
static TH2I * hDigiHit_RawPeakVsSlotID
jerror_t init(void)
Called once at program start.
static TH1I * hDigiHit_Pedestal
jerror_t fini(void)
Called after last event of last event source has been processed.
static TH1I * hDigiHit_PulseTime
static TH2F * hHit_TimeVsSlotID
static TH2F * hHit_HasTDCvsHasADC
static TH1I * hDigiHit_tdcadcTimeDiff
double beam_current
double getEhigh(unsigned int counter) const
static TH2I * hDigiHit_tdcTimeVsfadcTime
static int matches(char *d, char *c, t_iostream *fp)
Definition: hddm_t.c:356
const int Nslots
static TH1I * hHit_RawNHits
static TH1I * hDigiHit_tdcRawTime
static TH2I * hDigiHit_fadcTimeVsSlotID
static TH2I * hHit_tdcadcTimeDiffVsIntegral
static TH2I * hHit_HVidVsSlotID
static TH1I * hDigiHit_RawPeak
static TH1I * hDigiHit_NfadcHits_cut
static TProfile * hDigiHit_PedestalVsSlotID
static TH1I * hHit_NHits_ds
static TH1I * hDigiHit_fadcTime
static TH1I * hDigiHit_PulseNumber
static TH2I * hHit_tdcTimeVsSlotID
static TH2I * hDigiHit_NfadcHitsVsSlotID
static TH1I * hDigiHit_RawIntegral
static TH1I * hDigiHit_NfadcHits_multiPeak
static const unsigned int kCounterCount
Definition: DTAGHGeometry.h:28
static TH2I * hDigiHit_PulseNumberVsSlotID
const int NupstreamCounters
static TH2I * hDigiHit_tdcTimeVsSlotID
static TH2I * hDigiHit_fadcTimeVsQF_cut
static TH1I * hDigiHit_tdcTime
static TH1I * hDigiHit_QualityFactor
uint32_t pulse_number
pulse number for this channel, this event starting from 0
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
static TH1I * hHit_NHits
static TH2I * hHit_fadcTimeVsSlotID