Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JEventProcessor_TAGH_timewalk.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: JEventProcessor_TAGH_timewalk.cc
4 // Created: Fri Nov 13 10:13:02 EST 2015
5 // Creator: nsparks (on Linux cua2.jlab.org 3.10.0-229.20.1.el7.x86_64 x86_64)
6 //
7 
9 using namespace jana;
10 using namespace std;
11 
12 #include <TAGGER/DTAGHHit.h>
13 #include <RF/DRFTime.h>
14 
15 // Routine used to create our JEventProcessor
16 #include <JANA/JApplication.h>
17 #include <JANA/JFactory.h>
18 #include <TDirectory.h>
19 #include <TH1.h>
20 #include <TH2.h>
21 
22 const int Nslots = 274; // number of TAGH counter slots
26 
27 extern "C"{
28  void InitPlugin(JApplication *app){
29  InitJANAPlugin(app);
30  app->AddProcessor(new JEventProcessor_TAGH_timewalk());
31  }
32 } // "C"
33 
34 
35 //------------------
36 // JEventProcessor_TAGH_timewalk (Constructor)
37 //------------------
39 {
40 
41 }
42 
43 //------------------
44 // ~JEventProcessor_TAGH_timewalk (Destructor)
45 //------------------
47 {
48 
49 }
50 
51 //------------------
52 // init
53 //------------------
55 {
56  // This is called once at program startup. If you are creating
57  // and filling historgrams in this plugin, you should lock the
58  // ROOT mutex like this:
59 
60  const double Tl = -200.0;
61  const double Th = 600.0;
62  const int NTb = 8000;
63  TDirectory *mainDir = gDirectory;
64  TDirectory *taghDir = gDirectory->mkdir("TAGH_timewalk");
65  taghDir->cd();
66  gDirectory->mkdir("Offsets")->cd();
67  hTAGH_tdcadcTimeDiffVsSlotID = new TH2I("TAGH_tdcadcTimeDiffVsSlotID","TAGH TDC-ADC time difference vs. counter ID;counter ID;time(TDC) - time(ADC) [ns]",Nslots,0.5,0.5+Nslots,NTb,Tl,Th);
68  hTAGHRF_tdcTimeDiffVsSlotID = new TH2I("TAGHRF_tdcTimeDiffVsSlotID","TAGH-RF TDC time difference vs. counter ID;counter ID;time(TDC) - time(RF) [ns]",Nslots,0.5,0.5+Nslots,NTb,Tl,Th);
69  taghDir->cd();
70  gDirectory->mkdir("Timewalk")->cd();
71  for (int i = 0; i < 1+Nslots; i++) {
72  stringstream ss; ss << i;
73  TString name = "TAGHRF_tdcTimeDiffVsPulseHeight_" + ss.str();
74  TString title = "TAGH counter " + ss.str();
75  hTAGHRF_tdcTimeDiffVsPulseHeight[i] = new TH2I(name,title+";pulse height [fADC counts];time(TDC) - time(RF) [ns]",41,0,4100,50,-2.5,2.5);
76  }
77  mainDir->cd();
78 
79  return NOERROR;
80 }
81 
82 //------------------
83 // brun
84 //------------------
85 jerror_t JEventProcessor_TAGH_timewalk::brun(JEventLoop *eventLoop, int32_t runnumber)
86 {
87  // This is called whenever the run number changes
88  dRFTimeFactory = static_cast<DRFTime_factory*>(eventLoop->GetFactory("DRFTime"));
89  vector<const DRFTime*> locRFTimes;
90  eventLoop->Get(locRFTimes);
91  return NOERROR;
92 }
93 
94 //------------------
95 // evnt
96 //------------------
97 jerror_t JEventProcessor_TAGH_timewalk::evnt(JEventLoop *loop, uint64_t eventnumber)
98 {
99  // This is called for every event. Use of common resources like writing
100  // to a file or filling a histogram should be mutex protected. Using
101  // loop->Get(...) to get reconstructed objects (and thereby activating the
102  // reconstruction algorithm) should be done outside of any mutex lock
103  // since multiple threads may call this method at the same time.
104  vector<const DTAGHHit*> taghhits;
105  loop->Get(taghhits, "Calib");
106  const DRFTime* rfTime = nullptr;
107  vector <const DRFTime*> rfTimes;
108  loop->Get(rfTimes,"TAGH");
109  if (rfTimes.size() > 0)
110  rfTime = rfTimes[0];
111  else
112  return NOERROR;
113 
114  // FILL HISTOGRAMS
115  // Since we are filling histograms local to this plugin, it will not interfere with other ROOT operations: can use plugin-wide ROOT fill lock
116  japp->RootFillLock(this); //ACQUIRE ROOT FILL LOCK
117 
118  double t_RF = rfTime->dTime;
119  for (unsigned int i = 0; i < taghhits.size(); i++) {
120  const DTAGHHit *hit = taghhits[i];
121  if (!hit->has_TDC || !hit->has_fADC) continue;
122  int id = hit->counter_id;
123  double t_tdc = hit->t;
124  double t_adc = hit->time_fadc;
125  double pulse_height = hit->pulse_peak;
126  hTAGH_tdcadcTimeDiffVsSlotID->Fill(id,t_tdc-t_adc);
127  hTAGHRF_tdcTimeDiffVsSlotID->Fill(id,t_tdc-t_RF);
128  t_tdc = dRFTimeFactory->Step_TimeToNearInputTime(t_tdc,t_RF);
129  hTAGHRF_tdcTimeDiffVsPulseHeight[0]->Fill(pulse_height,t_tdc-t_RF);
130  hTAGHRF_tdcTimeDiffVsPulseHeight[id]->Fill(pulse_height,t_tdc-t_RF);
131  }
132 
133  japp->RootFillUnLock(this); //RELEASE ROOT FILL LOCK
134 
135  return NOERROR;
136 }
137 
138 //------------------
139 // erun
140 //------------------
142 {
143  // This is called whenever the run number changes, before it is
144  // changed to give you a chance to clean up before processing
145  // events from the next run number.
146  return NOERROR;
147 }
148 
149 //------------------
150 // fini
151 //------------------
153 {
154  // Called before program exit after event processing is finished.
155  return NOERROR;
156 }
static TH2I * hTAGHRF_tdcTimeDiffVsPulseHeight[1+Nslots]
double t
Definition: DTAGHHit.h:19
jerror_t fini(void)
Called after last event of last event source has been processed.
bool has_TDC
Definition: DTAGHHit.h:26
double dTime
Definition: DRFTime.h:24
TDirectory * mainDir
Definition: p2k_hists.C:2
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
jerror_t init(void)
Called once at program start.
JApplication * japp
int counter_id
Definition: DTAGHHit.h:20
InitPlugin_t InitPlugin
static TH2I * hTAGH_tdcadcTimeDiffVsSlotID
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
static TH2I * hTAGHRF_tdcTimeDiffVsSlotID
const int Nslots
bool has_fADC
Definition: DTAGHHit.h:26
double time_fadc
Definition: DTAGHHit.h:24
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
double pulse_peak
Definition: DTAGHHit.h:22
double Step_TimeToNearInputTime(double locTimeToStep, double locTimeToStepTo) const
DRFTime_factory * dRFTimeFactory