Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DEventProcessor_ccal_hits.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DEventProcessor_ccal_hits.cc
4 // Created: Mon Apr 3 11:38:03 EDT 2006
5 // Creator: davidl (on Darwin swire-b241.jlab.org 8.4.0 powerpc)
6 //
7 
8 #include <map>
9 using namespace std;
10 
12 
13 #include <DANA/DApplication.h>
14 #include <CCAL/DCCALDigiHit.h>
15 
16 #include <CCAL/DCCALHit.h>
17 
18 #include <DAQ/Df250WindowRawData.h>
19 #include <DAQ/Df250PulseData.h>
20 
21 
22 // Routine used to create our DEventProcessor
23 extern "C"{
24 void InitPlugin(JApplication *app){
25  InitJANAPlugin(app);
26  app->AddProcessor(new DEventProcessor_ccal_hits());
27 }
28 } // "C"
29 
30 
31 //------------------
32 // init
33 //------------------
35 {
36 
37  TDirectory *dir = new TDirectoryFile("CCAL","CCAL");
38  dir->cd();
39 
40  tree1 = new TTree( "digihit", "digihit" );
41 
42  tree1->Branch("nhit", &nhit, "nhit/I");
43  tree1->Branch("column", &column, "column[nhit]/I");
44  tree1->Branch("row", &row, "row[nhit]/I");
45  tree1->Branch("peak", &peak, "peak[nhit]/I");
46  tree1->Branch("integral", &integral, "integral[nhit]/I");
47  tree1->Branch("pedestal", &pedestal, "pedestal[nhit]/I");
48  tree1->Branch("time", &time, "time[nhit]/I");
49  tree1->Branch("qf", &qf, "qf[nhit]/I");
50 
51  tree1->Branch("waveform", &waveform, "waveform[nhit][100]/I");
52 
53  tree1->Branch("nsamples_integral", &nsamples_integral, "nsamples_integral/I");
54  tree1->Branch("nsamples_pedestal", &nsamples_pedestal, "nsamples_pedestal/I");
55 
56 
57 
58  for(Int_t ii = 0; ii < 12; ii++){
59  for(Int_t jj = 0; jj < 12; jj++){
60 
61  char title[30];
62 
63  int index = jj + ii*12;
64 
65  sprintf(title,"Waveform_%d",index);
66  ccal_wave[index] = new TProfile(title,title,100,-0.5,99.5,-10.,4096);
67 
68  sprintf(title,"Peak_%d",index);
69  ccal_peak[index] = new TH1F(title, title, 4096, -0.5, 4095.5);
70 
71  sprintf(title,"Int_%d",index);
72  ccal_int[index] = new TH1F(title, title, 500, 0., 20000.5);
73 
74  }
75  }
76 
77 
78  // Go back up to the parent directory
79  dir->cd("../");
80 
81 
82  return NOERROR;
83 }
84 
85 //------------------
86 // brun
87 //------------------
88 jerror_t DEventProcessor_ccal_hits::brun(JEventLoop *eventLoop, int32_t runnumber)
89 {
90 
91  return NOERROR;
92 }
93 
94 //------------------
95 // evnt
96 //------------------
97 jerror_t DEventProcessor_ccal_hits::evnt(JEventLoop *loop, uint64_t eventnumber)
98 {
99 
100  vector<const DCCALDigiHit*> ccal_digihits;
101 
102  vector<const DCCALHit*> ccal_hits;
103 
104  loop->Get(ccal_digihits);
105  loop->Get(ccal_hits);
106 
107  cout << " Event number = " << eventnumber << endl;
108  cout << " Number of digi hits = " << ccal_digihits.size() << endl;
109  cout << " Number of hits = " << ccal_hits.size() << endl;
110 
111 
112  nhit = 0;
113  memset(column,0,sizeof(column));
114  memset(row,0,sizeof(row));
115  memset(peak,0,sizeof(peak));
116  memset(integral,0,sizeof(integral));
117  memset(pedestal,0,sizeof(pedestal));
118  memset(time,0,sizeof(time));
119  memset(qf,0,sizeof(qf));
120  memset(waveform,0,sizeof(waveform));
121  nsamples_integral = 0;
122  nsamples_pedestal = 0;
123 
124 
125  japp->RootWriteLock();
126 
127 
128  for(unsigned int ii = 0; ii < ccal_digihits.size(); ii++){
129  const DCCALDigiHit *ccal_hit = ccal_digihits[ii];
130 
131  if(ii == 0){
132  nsamples_integral = ccal_hit->nsamples_integral;
133  nsamples_pedestal = ccal_hit->nsamples_pedestal;
134  }
135 
136  row[nhit] = ccal_hit->row;
137  column[nhit] = ccal_hit->column;
138 
139  peak[nhit] = ccal_hit->pulse_peak;
140  integral[nhit] = ccal_hit->pulse_integral;
141 
142  pedestal[nhit] = ccal_hit->pedestal;
143  time[nhit] = (ccal_hit->pulse_time & 0x7FC0) >> 6;
144  qf[nhit] = ccal_hit->QF;
145 
146  int index = ccal_hit->row*12 + ccal_hit->column;
147 
148  // Assume that baseline is at fadc count 100
149  ccal_peak[index]->Fill(float(ccal_hit->pulse_peak) - 100.);
150  ccal_int[index]->Fill(float(ccal_hit->pulse_integral) - 100*nsamples_integral);
151 
152 
153  const Df250WindowRawData *windorawdata;
154  const Df250PulseData *pulsedata;
155 
156  ccal_hit->GetSingle(pulsedata);
157 
158  if(pulsedata){;
159 
160  pulsedata->GetSingle(windorawdata);
161 
162  if(windorawdata){
163 
164  const vector<uint16_t> &samplesvector = windorawdata->samples;
165 
166  unsigned int nsamples = samplesvector.size();
167 
168  for(uint16_t samp = 0; samp < nsamples; samp++){
169  waveform[nhit][samp] = samplesvector[samp];
170 
171  ccal_wave[index]->Fill(float(samp),float(samplesvector[samp]));
172 
173  }
174 
175  }
176 
177  }
178 
179 
180  nhit++;
181 
182  }
183 
184 
185  if(nhit > 0){
186  tree1->Fill();
187  }
188 
189  japp->RootUnLock();
190 
191  return NOERROR;
192 }
193 
194 //------------------
195 // erun
196 //------------------
198 {
199  // Any final calculations on histograms (like dividing them)
200  // should be done here. This may get called more than once.
201  return NOERROR;
202 }
203 
204 //------------------
205 // fini
206 //------------------
208 {
209  return NOERROR;
210 }
211 
jerror_t evnt(JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
jerror_t brun(JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
uint32_t pulse_integral
identified pulse integral as returned by FPGA algorithm
Definition: DCCALDigiHit.h:20
sprintf(text,"Post KinFit Cut")
jerror_t init(void)
Called once at program start.
vector< uint16_t > samples
static char index(char c)
Definition: base64.cpp:115
JApplication * japp
uint32_t pedestal
pedestal info used by FPGA (if any)
Definition: DCCALDigiHit.h:22
InitPlugin_t InitPlugin
uint32_t nsamples_integral
number of samples used in integral
Definition: DCCALDigiHit.h:24
uint32_t nsamples_pedestal
number of samples used in pedestal
Definition: DCCALDigiHit.h:25
uint32_t QF
Quality Factor from FPGA algorithms.
Definition: DCCALDigiHit.h:23
static TH1I * pedestal[nChan]
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
uint32_t pulse_peak
maximum sample in pulse
Definition: DCCALDigiHit.h:26
TDirectory * dir
Definition: bcal_hist_eff.C:31
uint32_t pulse_time
identified pulse time as returned by FPGA algorithm
Definition: DCCALDigiHit.h:21
jerror_t fini(void)
Called after last event of last event source has been processed.