Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_PSC_TW.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_PSC_TW.cc
4 // Created: Fri Aug 21 10:42:28 EDT 2015
5 // Creator: aebarnes (on Linux ifarm1102 2.6.32-431.el6.x86_64 x86_64)
6 //
7 
9 using namespace jana;
10 
11 // Routine used to create our JEventProcessor
12 #include <JANA/JApplication.h>
13 #include <JANA/JFactory.h>
14 // ROOT header fiels
15 #include <TH2.h>
16 #include <TDirectory.h>
17 // RF header files
18 #include <RF/DRFTime_factory.h>
19 // PSC header files
21 
22 // Define constants
23 const uint32_t NMODULES = 8;
24 
25 const double TMIN = -5;
26 const double TMAX = 5;
27 const double TBINS = (TMAX - TMIN)/0.1;
28 
29 const double PMIN = 0;
30 const double PMAX = 2000;
31 const double PBINS = (PMAX - PMIN)/2;
32 
33 // Define histograms
34 static TH2F* h_dt_vs_pp_tdc_l[NMODULES];
35 static TH2F* h_dt_vs_pp_tdc_r[NMODULES];
36 static TH2F* h_dt_vs_pp_t_l[NMODULES];
37 static TH2F* h_dt_vs_pp_t_r[NMODULES];
38 
39 // Define variables
40 Int_t psc_mod_l;
41 Int_t psc_mod_r;
42 Double_t pp_l;
43 Double_t pp_r;
44 Double_t adc_l;
45 Double_t adc_r;
46 Double_t tdc_l;
47 Double_t tdc_r;
48 Double_t t_l;
49 Double_t t_r;
50 Double_t rf_l;
51 Double_t rf_r;
52 
53 // Define RFTime_factory
55 
56 extern "C"{
57 void InitPlugin(JApplication *app){
58  InitJANAPlugin(app);
59  app->AddProcessor(new JEventProcessor_PSC_TW());
60 }
61 } // "C"
62 
63 
64 //------------------
65 // JEventProcessor_PSC_TW (Constructor)
66 //------------------
68 {
69 
70 }
71 
72 //------------------
73 // ~JEventProcessor_PSC_TW (Destructor)
74 //------------------
76 {
77 
78 }
79 
80 //------------------
81 // init //------------------
83 {
84  // This is called once at program startup. If you are creating
85  // and filling historgrams in this plugin, you should lock the
86  // ROOT mutex like this:
87  //
88  // japp->RootWriteLock();
89  // ... fill historgrams or trees ...
90  // japp->RootUnLock();
91  //
92 
93  TDirectory *pscDir = gDirectory->mkdir("PSC_TW");
94  pscDir->cd();
95 
96  gDirectory->mkdir("tdc-rf")->cd();
97  for (uint32_t i = 0; i < NMODULES; ++i)
98  {
99  h_dt_vs_pp_tdc_l[i] = new TH2F(Form("h_dt_vs_pp_tdc_l_%i",i+1),
100  Form("#Delta t vs. pulse peak, raw TDC, left PSC %i;\
101  Pulse peak; #Delta t (raw TDC - RF)",i+1),
102  PBINS, PMIN, PMAX, TBINS, TMIN, TMAX);
103  h_dt_vs_pp_tdc_r[i] = new TH2F(Form("h_dt_vs_pp_tdc_r_%i",i+1),
104  Form("#Delta t vs. pulse peak, raw TDC, left PSC %i;\
105  Pulse peak; #Delta t (raw TDC - RF)",i+1),
106  PBINS, PMIN, PMAX, TBINS, TMIN, TMAX);
107  }
108  pscDir->cd();
109 
110  gDirectory->mkdir("t-rf")->cd();
111  for (uint32_t i = 0; i < NMODULES; ++i)
112  {
113  h_dt_vs_pp_t_l[i] = new TH2F(Form("h_dt_vs_pp_t_l_%i",i+1),
114  Form("#Delta t vs. pulse peak, corrected TDC, left PSC %i;\
115  Pulse Peak; #Delta t (raw TDC - RF)",i+1),
116  PBINS, PMIN, PMAX, TBINS, TMIN, TMAX);
117  h_dt_vs_pp_t_r[i] = new TH2F(Form("h_dt_vs_pp_t_r_%i",i+1),
118  Form("#Delta t vs. pulse peak, corrected TDC, left PSC %i;\
119  Pulse peak; #Delta t (raw TDC - RF)",i+1),
120  PBINS, PMIN, PMAX, TBINS, TMIN, TMAX);
121  }
122  pscDir->cd();
123 
124  return NOERROR;
125 }
126 
127 //------------------
128 // brun
129 //------------------
130 jerror_t JEventProcessor_PSC_TW::brun(JEventLoop *eventLoop, int32_t runnumber)
131 {
132  // This is called whenever the run number changes
133 
134  //////////////
135  // RF
136  //////////////
137 
138  return NOERROR;
139 }
140 
141 //------------------
142 // evnt
143 //------------------
144 jerror_t JEventProcessor_PSC_TW::evnt(JEventLoop *loop, uint64_t eventnumber)
145 {
146  // This is called for every event. Use of common resources like writing
147  // to a file or filling a histogram should be mutex protected. Using
148  // loop->Get(...) to get reconstructed objects (and thereby activating the
149  // reconstruction algorithm) should be done outside of any mutex lock
150  // since multiple threads may call this method at the same time.
151  // Here's an example:
152  //
153  // vector<const MyDataClass*> mydataclasses;
154  // loop->Get(mydataclasses);
155  //
156  // japp->RootWriteLock();
157  // ... fill historgrams or trees ...
158  // japp->RootUnLock();
159 
160  // Initialize RF time factory
161  auto locRFTimeFactory = static_cast<DRFTime_factory*>(loop->GetFactory("DRFTime"));
162 
163  vector<const DRFTime*> locRFTimes;
164  vector<const DPSCPair*> pairs;
165 
166  loop->Get(pairs);
167  loop->Get(locRFTimes,"PSC");
168  const DRFTime* locRFTime = NULL;
169 
170  if (locRFTimes.size() > 0)
171  locRFTime = locRFTimes[0];
172  else
173  return NOERROR;
174 
175  // Since we are filling histograms local to this plugin,
176  // it will not interfere with other ROOT operations:
177  // can use plugin-wide ROOT fill lock
178  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
179 
180  for (uint32_t i = 0; i < pairs.size(); ++i)
181  {
182  // Left modules
183  psc_mod_l = pairs[i]->ee.first->module;
184  pp_l = pairs[i]->ee.first->pulse_peak;
185  adc_l = pairs[i]->ee.first->time_fadc;
186  tdc_l = pairs[i]->ee.first->time_tdc;
187  t_l = pairs[i]->ee.first->t;
188  // Right modules
189  psc_mod_r = pairs[i]->ee.second->module;
190  pp_r = pairs[i]->ee.second->pulse_peak;
191  adc_r = pairs[i]->ee.second->time_fadc;
192  tdc_r = pairs[i]->ee.second->time_tdc;
193  t_r = pairs[i]->ee.second->t;
194 
195  // Use the ADC time to find the closest RF time.
198 
199  h_dt_vs_pp_tdc_l[psc_mod_l - 1]->Fill(pp_l, tdc_l - rf_l);
200  h_dt_vs_pp_t_l[psc_mod_l - 1]->Fill(pp_l, t_l - rf_l);
201  h_dt_vs_pp_tdc_r[psc_mod_r - 1]->Fill(pp_r, tdc_r - rf_r);
202  h_dt_vs_pp_t_r[psc_mod_r - 1]->Fill(pp_r, t_r - rf_r);
203  }
204 
205  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
206 
207  return NOERROR;
208 }
209 
210 //------------------
211 // erun
212 //------------------
214 {
215  // This is called whenever the run number changes, before it is
216  // changed to give you a chance to clean up before processing
217  // events from the next run number.
218  return NOERROR;
219 }
220 
221 //------------------
222 // fini
223 //------------------
225 {
226  // Called before program exit after event processing is finished.
227  return NOERROR;
228 }
Double_t tdc_l
Int_t psc_mod_r
const double PBINS
DRFTime_factory * locRFTimeFactory
Double_t adc_r
Double_t t_r
Double_t pp_l
Int_t psc_mod_l
Double_t pp_r
Double_t tdc_r
const double PMAX
static TH2F * h_dt_vs_pp_tdc_r[NMODULES]
double dTime
Definition: DRFTime.h:24
Double_t rf_l
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
static TH2F * h_dt_vs_pp_t_r[NMODULES]
JApplication * japp
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
static TH2F * h_dt_vs_pp_t_l[NMODULES]
const double PMIN
InitPlugin_t InitPlugin
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
static TH2F * h_dt_vs_pp_tdc_l[NMODULES]
const double TMAX
jerror_t init(void)
Called once at program start.
const uint32_t NMODULES
const double TMIN
Double_t rf_r
jerror_t fini(void)
Called after last event of last event source has been processed.
Double_t adc_l
const double TBINS
double Step_TimeToNearInputTime(double locTimeToStep, double locTimeToStepTo) const
Double_t t_l