Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_BCAL_saturation.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_BCAL_saturation.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 "BCAL/DBCALDigiHit.h"
18 #include "DANA/DStatusBits.h"
19 #include "DAQ/DEPICSvalue.h"
20 #include "DAQ/DF1TDCHit.h"
21 #include "DAQ/Df250PulseData.h"
22 #include "DAQ/Df250PulseIntegral.h"
23 #include "DAQ/Df250WindowRawData.h"
24 #include "TRIGGER/DL1Trigger.h"
25 
26 #include <TDirectory.h>
27 #include <TH2.h>
28 #include <TH1.h>
29 #include <TProfile2D.h>
30 #include <TStyle.h>
31 
32 
33 // root hist pointers
34 static TH1I *bcal_fadc_digi_integral = NULL;
35 static TH1I *bcal_fadc_digi_peak = NULL;
36 static TH2I *bcal_fadc_digi_integral_peak = NULL;
37 static TH1I *bcal_fadc_digi_QF = NULL;
38 static TH1I *bcal_fadc_digi_time = NULL;
39 static TH1I *bcal_fadc_digi_time_saturate = NULL;
40 static TH1I *bcal_saturated_counter = NULL;
41 
42 static TH2F *bcalUS_waveform[2][4] = {NULL,NULL};
43 static TH2F *bcalDS_waveform[2][4] = {NULL,NULL};
44 
45 static TH1I *bcal_num_events;
46 
47 //----------------------------------------------------------------------------------
48 
49 
50 // Routine used to create our JEventProcessor
51 extern "C"{
52  void InitPlugin(JApplication *app){
53  InitJANAPlugin(app);
54  app->AddProcessor(new JEventProcessor_BCAL_saturation());
55  }
56 }
57 
58 
59 //----------------------------------------------------------------------------------
60 
61 
63 }
64 
65 
66 //----------------------------------------------------------------------------------
67 
68 
70 }
71 
72 
73 //----------------------------------------------------------------------------------
74 
76 
77  // First thread to get here makes all histograms. If one pointer is
78  // already not NULL, assume all histograms are defined and return now
79  if(bcal_fadc_digi_integral != NULL){
80  return NOERROR;
81  }
82 
83  // create root folder for bcal and cd to it, store main dir
84  TDirectory *main = gDirectory;
85  gDirectory->mkdir("bcal_saturation")->cd();
86 
87  gStyle->SetTitleOffset(1, "Y");
88  gStyle->SetTitleSize(0.05,"xyz");
89  gStyle->SetTitleSize(0.08,"h");
90  gStyle->SetLabelSize(0.05,"xyz");
91  gStyle->SetTitleX(0);
92  gStyle->SetTitleAlign(13);
93  gStyle->SetNdivisions(505,"xy");
94 
95  bcal_num_events = new TH1I("bcal_num_events","BCAL Number of events",1, 0.5, 1.5);
96 
97  bcal_saturated_counter = new TH1I("bcal_saturated_counter", "BCAL # saturated fADCs per event; # saturated fADCs", 10, 0, 10);
98 
99  TString locSaturateName[2] = {"saturate", "noSaturate"};
100  for(int locSaturate = 0; locSaturate<2; locSaturate++){
101  for(int locLayer=0; locLayer<4; locLayer++){
102 
103  bcalUS_waveform[locSaturate][locLayer] = new TH2F(Form("bcalUS_waveform_%s_layer%d", locSaturateName[locSaturate].Data(), locLayer+1), "Upstream BCAL Waveforms; DigiHit Number; Samples", 10000, 0, 10000, 100, 0, 100);
104  bcalDS_waveform[locSaturate][locLayer] = new TH2F(Form("bcalDS_waveform_%s_layer%d", locSaturateName[locSaturate].Data(), locLayer+1), "Downstream BCAL Waveforms; DigiHit Number; Samples", 10000, 0, 10000, 100, 0, 100);
105  }
106  }
107 
108  bcal_fadc_digi_integral = new TH1I("bcal_fadc_digi_integral","BCAL Integral (DBCALDigiHit);Integral (fADC counts)", 600, 0, 120000);
109  bcal_fadc_digi_peak = new TH1I("bcal_fadc_digi_peak","BCAL Peak (DBCALDigiHit);Peak (fADC counts)", 1000, 0, 5000);
110  bcal_fadc_digi_integral_peak = new TH2I("bcal_fadc_digi_integral_peak","BCAL Integral vs Peak (DBCALDigiHit);Peak (fADC counts);Integral (fADC counts)", 1000, 0, 5000, 600, 0, 120000);
111  bcal_fadc_digi_QF = new TH1I("bcal_fadc_digi_QF","Qualtiy Factor (DBCALDigiHit);Qualtiy Factor", 20, 0, 20);
112  bcal_fadc_digi_time = new TH1I("bcal_fadc_digi_time","ADC Time (DBCALDigiHit);Time (fADC time/62.5 ps)", 550, -600, 6000);
113  bcal_fadc_digi_time_saturate = new TH1I("bcal_fadc_digi_time_saturate","Saturated pulse ADC Time (DBCALDigiHit);Time (fADC time/62.5 ps)", 550, -600, 6000);
114 
115  // back to main dir
116  main->cd();
117 
118  for(int locSaturate = 0; locSaturate<2; locSaturate++){
119  for(int locLayer=0; locLayer<4; locLayer++){
120  waveformCounterUS[locSaturate][locLayer] = 0;
121  waveformCounterDS[locSaturate][locLayer] = 0;
122  }
123  }
124 
125  return NOERROR;
126 }
127 
128 
129 //----------------------------------------------------------------------------------
130 
131 
132 jerror_t JEventProcessor_BCAL_saturation::brun(JEventLoop *eventLoop, int32_t runnumber) {
133  // This is called whenever the run number changes
134  return NOERROR;
135 }
136 
137 
138 //----------------------------------------------------------------------------------
139 
140 
141 jerror_t JEventProcessor_BCAL_saturation::evnt(JEventLoop *loop, uint64_t eventnumber) {
142  // This is called for every event. Use of common resources like writing
143  // to a file or filling a histogram should be mutex protected. Using
144  // loop-Get(...) to get reconstructed objects (and thereby activating the
145  // reconstruction algorithm) should be done outside of any mutex lock
146  // since multiple threads may call this method at the same time.
147 
148  vector<const DBCALDigiHit*> dbcaldigihits;
149 
150  // First check that this is not a font panel trigger or no trigger
151  bool goodtrigger=1;
152 
153  const DL1Trigger *trig = NULL;
154  try {
155  loop->GetSingle(trig);
156  } catch (...) {}
157  if (trig) {
158  if (trig->fp_trig_mask){
159  goodtrigger=0;
160  }
161  } else {
162  // HDDM files are from simulation, so keep them even though they have no trigger
163  bool locIsHDDMEvent = loop->GetJEvent().GetStatusBit(kSTATUS_HDDM);
164  if (!locIsHDDMEvent) goodtrigger=0;
165  }
166 
167  if (!goodtrigger) {
168  return NOERROR;
169  }
170 
171  loop->Get(dbcaldigihits);
172 
173  int locSaturatedCounter = 0;
174 
175  // FILL HISTOGRAMS
176  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
177  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
178 
179  if( (dbcaldigihits.size() > 0) )
180  bcal_num_events->Fill(1);
181 
182  // Digitized fADC hits for bcal
183  for(unsigned int i=0; i<dbcaldigihits.size(); i++) {
184  const DBCALDigiHit *hit = dbcaldigihits[i];
185  int locLayer = hit->layer-1;
186 
188  if(hit->pulse_peak) {
189  bcal_fadc_digi_peak->Fill(hit->pulse_peak);
191  }
192  vector<const Df250WindowRawData*> f250WindowRawData;
193  hit->Get(f250WindowRawData);
194  if(hit->pulse_peak > 3000 && f250WindowRawData.size() > 0) {
195  // Get a vector of the samples for this channel
196  const vector<uint16_t> &samples = f250WindowRawData[0]->samples;
197  uint nsamples=samples.size();
198  for(uint isample=0; isample<nsamples; isample++) {
199  if(hit->pulse_peak > 4094)
200  if(hit->end == DBCALGeometry::kUpstream)
201  bcalUS_waveform[0][locLayer]->Fill(waveformCounterUS[0][locLayer], isample, samples[isample]);
202  else
203  bcalDS_waveform[0][locLayer]->Fill(waveformCounterDS[0][locLayer], isample, samples[isample]);
204  else
205  if(hit->end == DBCALGeometry::kUpstream)
206  bcalUS_waveform[1][locLayer]->Fill(waveformCounterUS[1][locLayer], isample, samples[isample]);
207  else
208  bcalDS_waveform[1][locLayer]->Fill(waveformCounterDS[1][locLayer], isample, samples[isample]);
209  }
210  if(hit->pulse_peak > 4094)
211  if(hit->end == DBCALGeometry::kUpstream)
212  waveformCounterUS[0][locLayer]++;
213  else
214  waveformCounterDS[0][locLayer]++;
215  else
216  if(hit->end == DBCALGeometry::kUpstream)
217  waveformCounterUS[1][locLayer]++;
218  else
219  waveformCounterDS[1][locLayer]++;
220  }
221 
222  bcal_fadc_digi_time->Fill(hit->pulse_time);
223  if(hit->pulse_peak > 4094) {
225  locSaturatedCounter++;
226  }
227  }
228  bcal_saturated_counter->Fill(locSaturatedCounter);
229 
230  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
231 
232  return NOERROR;
233 }
234 
235 
236 //----------------------------------------------------------------------------------
237 
238 
240  // This is called whenever the run number changes, before it is
241  // changed to give you a chance to clean up before processing
242  // events from the next run number.
243 
244  return NOERROR;
245 }
246 
247 
248 //----------------------------------------------------------------------------------
249 
250 
252  // Called before program exit after event processing is finished.
253  return NOERROR;
254 }
255 
256 
257 //----------------------------------------------------------------------------------
258 //----------------------------------------------------------------------------------
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
uint32_t fp_trig_mask
Definition: DL1Trigger.h:19
static TH2F * bcalDS_waveform[2][4]
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
static TH1I * bcal_saturated_counter
static TH1I * bcal_fadc_digi_time
trig[33-1]
JApplication * japp
DBCALGeometry::End end
Definition: DBCALDigiHit.h:28
uint32_t pulse_peak
identified pulse height as returned by FPGA algorithm
Definition: DBCALDigiHit.h:30
static TH1I * bcal_fadc_digi_QF
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
uint32_t pulse_integral
identified pulse integral as returned by FPGA algorithm
Definition: DBCALDigiHit.h:29
jerror_t init(void)
Called once at program start.
InitPlugin_t InitPlugin
static TH2I * bcal_fadc_digi_integral_peak
static TH2F * bcalUS_waveform[2][4]
static TH1I * bcal_num_events
static TH1I * bcal_fadc_digi_time_saturate
uint32_t pulse_time
identified pulse time as returned by FPGA algorithm
Definition: DBCALDigiHit.h:31
static TH1I * bcal_fadc_digi_peak
jerror_t fini(void)
Called after last event of last event source has been processed.
static TH1I * bcal_fadc_digi_integral
int main(int argc, char *argv[])
Definition: gendoc.cc:6