Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_st_tw_corr_auto.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_st_tw_corr_auto.cc
4 // Created: Mon Oct 26 10:35:45 EDT 2015
5 // Creator: mkamel (on Linux ifarm1401 2.6.32-431.el6.x86_64 x86_64)
6 //
8 #include "TRIGGER/DTrigger.h"
9 
10 // Routine used to create our JEventProcessor
11 #include <JANA/JApplication.h>
12 #include <JANA/JFactory.h>
13 
14 extern "C"{
15 void InitPlugin(JApplication *app){
16  InitJANAPlugin(app);
17  app->AddProcessor(new JEventProcessor_st_tw_corr_auto());
18 }
19 } // "C"
20 
21 //------------------
22 // JEventProcessor_st_tw_corr_auto (Constructor)
23 //------------------
25 {
26 
27 }
28 
29 //------------------
30 // ~JEventProcessor_st_tw_corr_auto (Destructor)
31 //------------------
33 {
34 
35 }
36 
37 //------------------
38 // init
39 //------------------
41 {
42  // This is called once at program startup. If you are creating
43  // and filling historgrams in this plugin, you should lock the
44  // ROOT mutex like this:
45  //
46  // japp->RootWriteLock();
47  // ... fill historgrams or trees ...
48  // japp->RootUnLock();
49  //
51  gPARMS->SetDefaultParameter("SC:USE_TIMEWALK_CORRECTION", USE_TIMEWALK_CORRECTION,
52  "Flag to decide if timewalk corrections should be applied.");
53  // ***************** define constants*************************
54  NCHANNELS = 30;
55  tdc_thresh_mV = 50.0;
56  tdc_gain_factor = 5.0;
57  adc_max_chan = 4096.0;
58  adc_max_mV = 2000.0;
60  // **************** define histograms *************************
61 
62  h_pp_chan = new TH1I*[NCHANNELS];
63  h_stt_chan = new TH1I*[NCHANNELS];
64  h2_stt_vs_pp_chan = new TH2I*[NCHANNELS];
65  h1_st_corr_time = new TH1I*[NCHANNELS];
66  h2_st_corr_vs_pp = new TH2I*[NCHANNELS];
67  // Book Histos
68  for (Int_t i = 0; i < NCHANNELS; i++)
69  {
70  h2_stt_vs_pp_chan[i] = new TH2I(Form("stt_vs_pp_chan_%i", i+1), "Hit Time vs. Pulse Peak; Pulse Peak (channels); #delta_{t} (ns)", 4096,0.0,4095.0, 160, -5.0, 5.0);
71 
72  h_pp_chan[i] = new TH1I(Form("pp_chan_%i", i+1), "Pulse Peak; Pulse Peak (channels); Counts", 4096,0.0,4095.0);
73 
74  h_stt_chan[i] = new TH1I(Form("h_stt_chan_%i",i+1), "ST Time; ST Time (ns); Counts", 160, -5.0, 5.0);
75 
76  h1_st_corr_time[i] = new TH1I(Form("h1_st_corr_time_%i",i+1), "ST Time; ST Time (ns); Counts", 160, -5.0, 5.0);
77 
78  h2_st_corr_vs_pp[i]= new TH2I(Form("h2_st_corr_vs_pp_%i", i+1), "Hit Time vs. Pulse Peak; Pulse Peak (channels); #delta_{t} (ns)", 4096,0.0,4095.0, 160, -5.0, 5.0);
79  }
80 
81  return NOERROR;
82 }
83 //------------------
84 // brun
85 //------------------
86 jerror_t JEventProcessor_st_tw_corr_auto::brun(JEventLoop *eventLoop, int32_t runnumber)
87 {
88  // This is called whenever the run number changes
89  // load constant tables
90  // timewalk_parameters (timewalk_parms)
91  if(eventLoop->GetCalib("START_COUNTER/timewalk_parms_v2", timewalk_parameters))
92  jout << "Error loading /START_COUNTER/timewalk_parms_v2 !" << endl;
93 
94  return NOERROR;
95 }
96 
97 //------------------
98 // evnt
99 //------------------
100 jerror_t JEventProcessor_st_tw_corr_auto::evnt(JEventLoop *loop, uint64_t eventnumber)
101 {
102  // This is called for every event. Use of common resources like writing
103  // to a file or filling a histogram should be mutex protected. Using
104  // loop->Get(...) to get reconstructed objects (and thereby activating the
105  // reconstruction algorithm) should be done outside of any mutex lock
106  // since multiple threads may call this method at the same time.
107  // Here's an example:
108  //
109  // vector<const MyDataClass*> mydataclasses;
110  // loop->Get(mydataclasses);
111  //
112  // japp->RootWriteLock();
113  // ... fill historgrams or trees ...
114  // japp->RootUnLock();
115 
116  // select events with physics events, i.e., not LED and other front panel triggers
117  const DTrigger* locTrigger = NULL;
118  loop->GetSingle(locTrigger);
119  if(locTrigger->Get_L1FrontPanelTriggerBits() != 0)
120  return NOERROR;
121 
122  vector<const DSCHit*> st_hits;
123  loop->Get(st_hits);
124 
125  for (unsigned int k = 0; k < st_hits.size(); k++)
126  {
127  if(!st_hits[k]->has_fADC || !st_hits[k]->has_TDC)
128  continue;
129  // get the associated digihit so that we can access the
130  // pulse peak and pedestal information
131  const DSCDigiHit* the_digihit = NULL;
132  st_hits[k]->GetSingle(the_digihit);
133  if(the_digihit == NULL)
134  continue;
135 
136  // Extract the information we're interested in
137  double pulse_peak = the_digihit->pulse_peak;
138  double pedestal = the_digihit->pedestal;
139  double nsamples_pedestal = the_digihit->nsamples_pedestal;
140  double adc_pp = pulse_peak - pedestal/nsamples_pedestal;
141 
142  double T = st_hits[k]->t; // this is always set to the TDC time
143  // timewalk corrections are controlled by command line flag
144  double adc_t = st_hits[k]->t_fADC;
145  int sector = st_hits[k]->sector;
146 
147  // FILL HISTOGRAMS
148  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
149  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
150  {
151  double st_time = T - adc_t;
152  h_stt_chan[sector-1]->Fill(st_time);
154  {
155  h_pp_chan[sector-1]->Fill(adc_pp);
156  h2_stt_vs_pp_chan[sector-1]->Fill(adc_pp,st_time);
157  }
158  // Correct for the time walk
159  else if (USE_TIMEWALK_CORRECTION==1)
160  {
161  double C0 = timewalk_parameters[sector -1][0];
162  double C1 = timewalk_parameters[sector -1][1];
163  double C2 = timewalk_parameters[sector -1][2];
164  double A_THRESH = timewalk_parameters[sector -1][3];
165  double A0 = timewalk_parameters[sector -1][4];
166  // cout<< C1 << "\t" << C2 << "\t" << A_THRESH << "\t" << A0 <<endl;
167 
168  st_time -= C0 + C1 *(TMath::Power(A0/A_THRESH, C2));
169  // cout << "st_time =" << st_time << endl;
170  h1_st_corr_time[sector-1]->Fill(st_time);
171  h2_st_corr_vs_pp[sector-1]->Fill(adc_pp,st_time);
172  }
173  }
174  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
175  }
176  return NOERROR;
177 }
178 
179 //------------------
180 // erun
181 //------------------
183 {
184  // This is called whenever the run number changes, before it is
185  // changed to give you a chance to clean up before processing
186  // events from the next run number.
187  return NOERROR;
188 }
189 
190 //------------------
191 // fini
192 //------------------
194 {
195  return NOERROR;
196 }
197 
uint32_t pedestal
pedestal info used by FPGA (if any)
Definition: DSCDigiHit.h:21
#define A0
Definition: src/md5.cpp:22
#define C0
Definition: src/md5.cpp:24
uint32_t Get_L1FrontPanelTriggerBits(void) const
JApplication * japp
jerror_t fini(void)
Called after last event of last event source has been processed.
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
InitPlugin_t InitPlugin
uint32_t nsamples_pedestal
number of samples used in pedestal
Definition: DSCDigiHit.h:24
jerror_t init(void)
Called once at program start.
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
static TH1I * pedestal[nChan]
uint32_t pulse_peak
maximum sample in pulse
Definition: DSCDigiHit.h:25
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.