Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_TAGM_TW.cc
Go to the documentation of this file.
1 // Plugin for the tagger microscope time-walk corrections
2 // Author: aebarnes
3 
5 using namespace jana;
6 
7 // Routine used to create our JEventProcessor
8 #include <JANA/JApplication.h>
9 #include <JANA/JFactory.h>
10 // ROOT header files
11 #include <TH1.h>
12 #include <TH2.h>
13 #include <TDirectory.h>
14 // GlueX header files
15 #include "TAGGER/DTAGMHit.h"
16 #include "RF/DRFTime_factory.h"
17 
18 // Define constants
19 const uint32_t NCOLUMNS = 102;
20 const uint32_t NSINGLES = 20;
21 const uint32_t NROWS = 5;
22 
23 const int32_t PMIN = 0;
24 const int32_t PMAX = 2000;
25 const int32_t PBIN = (PMAX - PMIN)/16;
26 
27 const int32_t TMIN = -50;
28 const int32_t TMAX = 50;
29 const int32_t TBIN = (TMAX - TMIN)/0.05;
30 
31 const double TMIN_RF = -2.0;
32 const double TMAX_RF = 2.0;
33 const double TBIN_RF = (TMAX_RF - TMIN_RF)/0.05;
34 
35 const double TMIN_TW = -10.0;
36 const double TMAX_TW = 15.0;
37 const double TBIN_TW = (TMAX_TW - TMIN_TW)/0.05;
38 
39 // Define histograms
40 // timewalk
41 static TH2I* h_dt_vs_pp[NCOLUMNS];
42 static TH2I* h_dt_vs_pp_tdc[NCOLUMNS];
43 static TH2I* h_dt_vs_pp_ind[NROWS][4];
44 static TH2I* h_dt_vs_pp_tdc_ind[NROWS][4];
45 // tdc - adc
46 static TH2I* h_tdc_adc_all;
47 static TH2I* h_tdc_adc_all_ind;
48 static TH2I* h_t_adc_all;
49 static TH2I* h_t_adc_all_ind;
50 // adc - rf
51 static TH2I* h_adc_rf_all;
52 static TH2I* h_adc_rf_all_ind;
53 
54 // Define RFTime_factory
56 
57 extern "C"{
58 void InitPlugin(JApplication *app){
59  InitJANAPlugin(app);
60  app->AddProcessor(new JEventProcessor_TAGM_TW());
61 }
62 } // "C"
63 
64 
65 //------------------
66 // JEventProcessor_TAGM_TW (Constructor)
67 //------------------
69 {
70 
71 }
72 
73 //------------------
74 // ~JEventProcessor_TAGM_TW (Destructor)
75 //------------------
77 {
78 
79 }
80 
81 //------------------
82 // init
83 //------------------
85 {
86  TDirectory *main = gDirectory;
87 
88  TDirectory *tagmDir = gDirectory->mkdir("TAGM_TW");
89  tagmDir->cd();
90 
91  gDirectory->mkdir("tdc-rf")->cd();
92  for (uint32_t i = 0; i < NCOLUMNS; ++i)
93  {
94  h_dt_vs_pp_tdc[i] = new TH2I(Form("h_dt_vs_pp_tdc_%i", i+1),
95  Form("#delta t vs. pulse peak for TAGM column %i;\
96  Pulse peak (adc counts);TDC - RF (ns)", i+1), \
97  PBIN, PMIN, PMAX, TBIN, TMIN, TMAX);
98  }
99  for (uint32_t i = 0; i < NROWS; ++i)
100  {
101  for (uint32_t j = 0; j < 4; ++j)
102  {
103  h_dt_vs_pp_tdc_ind[i][j] = new TH2I(Form("h_dt_vs_pp_tdc_ind_%i_%i", i+1, j+1),
104  Form("#delta t vs. pulse peak for TAGM ind. column %i, row %i;\
105  Pulse peak (adc counts);TDC - RF (ns)", j+1, i+1),
106  PBIN, PMIN, PMAX, TBIN, TMIN, TMAX);
107 
108  }
109  }
110  tagmDir->cd();
111 
112  gDirectory->mkdir("t-rf")->cd();
113  for (uint32_t i = 0; i < NCOLUMNS; ++i)
114  {
115  h_dt_vs_pp[i] = new TH2I(Form("h_dt_vs_pp_%i", i+1),
116  Form("#delta t vs. pulse peak for TAGM column %i;\
117  Pulse peak (adc counts);T - RF (ns)", i+1),\
119  }
120  for (uint32_t i = 0; i < NROWS; ++i)
121  {
122  for (uint32_t j = 0; j < 4; ++j)
123  {
124  h_dt_vs_pp_ind[i][j] = new TH2I(Form("h_dt_vs_pp_ind_%i_%i", i+1, j+1),
125  Form("#delta t vs. pulse peak for TAGM ind. column %i, row %i;\
126  Pulse peak (adc counts);T - RF (ns)", j+1, i+1),
128 
129  }
130  }
131  tagmDir->cd();
132 
133  h_tdc_adc_all = new TH2I("tdc_adc_all","Summed channels;TDC - ADC (ns);Column",
134  TBIN/10, TMIN, TMAX, NCOLUMNS,1, NCOLUMNS+1);
135  h_tdc_adc_all_ind = new TH2I("tdc_adc_all_ind","Individual channels;TDC - ADC (ns);Column",
136  TBIN/10, TMIN, TMAX, NSINGLES,1, NSINGLES+1);
137  h_t_adc_all = new TH2I("t_adc_all","Summed channels;TDC - ADC (ns);Column",
138  TBIN/10, TMIN, TMAX, NCOLUMNS,1, NCOLUMNS+1);
139  h_t_adc_all_ind = new TH2I("t_adc_all_ind","Individual channels;TDC - ADC (ns);Column",
140  TBIN/10, TMIN, TMAX, NSINGLES,1, NSINGLES+1);
141  h_adc_rf_all = new TH2I("adc_rf_all","Summed channels;ADC - RF (ns);Column",
142  TBIN_RF, TMIN_RF, TMAX_RF, NCOLUMNS,1, NCOLUMNS+1);
143  h_adc_rf_all_ind = new TH2I("adc_rf_all_ind","Individual channels;ADC - RF (ns);Column",
145 
146  main->cd();
147 
148  return NOERROR;
149 
150 }
151 
152 //------------------
153 // brun
154 //------------------
155 jerror_t JEventProcessor_TAGM_TW::brun(JEventLoop *eventLoop, int32_t runnumber)
156 {
157  // Initialize RF time factory
158  dRFTimeFactory = static_cast<DRFTime_factory*>(eventLoop->GetFactory("DRFTime"));
159 
160  // be sure that DRFTime_factory::init() and brun() are called
161  vector<const DRFTime*> locRFTimes;
162  eventLoop->Get(locRFTimes);
163 
164  return NOERROR;
165 
166 }
167 
168 //------------------
169 // evnt
170 //------------------
171 jerror_t JEventProcessor_TAGM_TW::evnt(JEventLoop *loop, uint64_t eventnumber)
172 {
173  // FILL HISTOGRAMS
174  // Since we are filling histograms local to this plugin, it will not interfere
175  // with other ROOT operations: can use plugin-wide ROOT fill lock
176 
177  vector<const DTAGMHit*> hits;
178  loop->Get(hits, "Calib");
179 
180  vector<const DRFTime*> locRFTimes;
181  loop->Get(locRFTimes, "TAGH");
182  const DRFTime* locRFTime = NULL;
183 
184  if (locRFTimes.size() > 0)
185  locRFTime = locRFTimes[0];
186  else
187  return NOERROR;
188 
189  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
190 
191  for (uint32_t i = 0; i < hits.size(); ++i) {
192  if (!hits[i]->has_TDC || !hits[i]->has_fADC) continue;
193  int col = hits[i]->column;
194  int row = hits[i]->row;
195  double pp = hits[i]->pulse_peak;
196  double tm_t = hits[i]->t;
197  double adc_t = hits[i]->time_fadc;
198  double tdc_t = hits[i]->time_tdc;
199  double rf_adc = dRFTimeFactory->Step_TimeToNearInputTime(locRFTime->dTime, adc_t);
200 
201  if (row == 0)
202  {
203  h_dt_vs_pp[col-1]->Fill(pp, tm_t - rf_adc);
204  h_dt_vs_pp_tdc[col-1]->Fill(pp, tdc_t - rf_adc);
205  h_tdc_adc_all->Fill(tdc_t - adc_t, col);
206  h_t_adc_all->Fill(tm_t - adc_t, col);
207  h_adc_rf_all->Fill(adc_t - rf_adc, col);
208  }
209  else
210  {
211  if (col == 9)
212  {
213  h_dt_vs_pp_ind[row-1][0]->Fill(pp,tm_t - rf_adc);
214  h_dt_vs_pp_tdc_ind[row-1][0]->Fill(pp,tdc_t - rf_adc);
215  h_tdc_adc_all_ind->Fill(tdc_t - adc_t, row);
216  h_t_adc_all_ind->Fill(tm_t - adc_t, row);
217  h_adc_rf_all_ind->Fill(adc_t - rf_adc, row);
218  }
219  else if (col == 27)
220  {
221  h_dt_vs_pp_ind[row-1][1]->Fill(pp,tm_t - rf_adc);
222  h_dt_vs_pp_tdc_ind[row-1][1]->Fill(pp,tdc_t - rf_adc);
223  h_tdc_adc_all_ind->Fill(tdc_t - adc_t, row + 5);
224  h_t_adc_all_ind->Fill(tm_t - adc_t, row + 5);
225  h_adc_rf_all_ind->Fill(adc_t - rf_adc, row + 5);
226  }
227  else if (col == 81)
228  {
229  h_dt_vs_pp_ind[row-1][2]->Fill(pp,tm_t - rf_adc);
230  h_dt_vs_pp_tdc_ind[row-1][2]->Fill(pp,tdc_t - rf_adc);
231  h_tdc_adc_all_ind->Fill(tdc_t - adc_t, row + 10);
232  h_t_adc_all_ind->Fill(tm_t - adc_t, row + 10);
233  h_adc_rf_all_ind->Fill(adc_t - rf_adc, row + 10);
234  }
235  else if (col == 99)
236  {
237  h_dt_vs_pp_ind[row-1][3]->Fill(pp,tm_t - rf_adc);
238  h_dt_vs_pp_tdc_ind[row-1][3]->Fill(pp,tdc_t - rf_adc);
239  h_tdc_adc_all_ind->Fill(tdc_t - adc_t, row + 15);
240  h_t_adc_all_ind->Fill(tm_t - adc_t, row + 15);
241  h_adc_rf_all_ind->Fill(adc_t - rf_adc, row + 15);
242  }
243  }
244  }
245 
246  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
247 
248  return NOERROR;
249 }
250 
251 //------------------
252 // erun
253 //------------------
255 {
256  // This is called whenever the run number changes, before it is
257  // changed to give you a chance to clean up before processing
258  // events from the next run number.
259  return NOERROR;
260 }
261 
262 //------------------
263 // fini
264 //------------------
266 {
267  // Called before program exit after event processing is finished.
268  return NOERROR;
269 }
const uint32_t NROWS
jerror_t init(void)
Called once at program start.
jerror_t fini(void)
Called after last event of last event source has been processed.
const int32_t PBIN
static TH2I * h_dt_vs_pp_tdc_ind[NROWS][4]
const double PMAX
static TH2I * h_tdc_adc_all
static TH2I * h_dt_vs_pp_tdc[NCOLUMNS]
double dTime
Definition: DRFTime.h:24
const double TBIN_TW
const double TMIN_TW
JApplication * japp
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
TH2I * h_dt_vs_pp
Definition: psc_mon.C:14
const double PMIN
InitPlugin_t InitPlugin
static TH2I * h_dt_vs_pp_ind[NROWS][4]
const double TBIN_RF
const double TMAX_RF
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
static TH2I * h_tdc_adc_all_ind
const uint32_t NCOLUMNS
const double TMAX_TW
const double TMAX
const int32_t TBIN
static TH2I * h_adc_rf_all
const double TMIN_RF
const double TMIN
static TH2I * h_t_adc_all
double Step_TimeToNearInputTime(double locTimeToStep, double locTimeToStepTo) const
const uint32_t NSINGLES
DRFTime_factory * dRFTimeFactory
static TH2I * h_t_adc_all_ind
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
int main(int argc, char *argv[])
Definition: gendoc.cc:6
static TH2I * h_adc_rf_all_ind