Hall C ROOT/C++ Analyzer (hcana)
THcScintillatorPlane.cxx
Go to the documentation of this file.
1 
9 #include "TMath.h"
10 #include "THcScintillatorPlane.h"
11 #include "THcScintPlaneCluster.h"
12 #include "TClonesArray.h"
13 #include "THcSignalHit.h"
14 #include "THcHodoHit.h"
15 #include "THcGlobals.h"
16 #include "THcParmList.h"
17 #include "THcHitList.h"
18 #include "THcHodoscope.h"
19 #include "TClass.h"
20 #include "THcRawAdcHit.h"
21 #include "THcRawTdcHit.h"
22 
23 #include <cstring>
24 #include <cstdio>
25 #include <cstdlib>
26 #include <iostream>
27 
28 using namespace std;
29 
31 
32 //______________________________________________________________________________
34  const char* description,
35  const Int_t planenum,
36  THaDetectorBase* parent )
37 : THaSubDetector(name,description,parent),
38  fParentHitList(0), fCluster(0),frPosAdcErrorFlag(0), frNegAdcErrorFlag(0),
39  frPosTDCHits(0), frNegTDCHits(0), frPosADCHits(0), frNegADCHits(0),
40  frPosADCSums(0), frNegADCSums(0), frPosADCPeds(0), frNegADCPeds(0),
41  fHodoHits(0), frPosTdcTimeRaw(0), frPosAdcPedRaw(0),
42  frPosAdcPulseIntRaw(0), frPosAdcPulseAmpRaw(0),
43  frPosAdcPulseTimeRaw(0), frPosTdcTime(0), frPosAdcPed(0),
44  frPosAdcPulseInt(0), frPosAdcPulseAmp(0), frPosAdcPulseTime(0),
45  frNegTdcTimeRaw(0), frNegAdcPedRaw(0), frNegAdcPulseIntRaw(0),
46  frNegAdcPulseAmpRaw(0), frNegAdcPulseTimeRaw(0), frNegTdcTime(0),
47  frNegAdcPed(0), frNegAdcPulseInt(0), frNegAdcPulseAmp(0),
48  frNegAdcPulseTime(0), fPosCenter(0), fHodoPosMinPh(0),
49  fHodoNegMinPh(0), fHodoPosPhcCoeff(0), fHodoNegPhcCoeff(0),
50  fHodoPosTimeOffset(0), fHodoNegTimeOffset(0), fHodoVelLight(0),
51  fHodoPosInvAdcOffset(0), fHodoNegInvAdcOffset(0),
52  fHodoPosAdcTimeWindowMin(0), fHodoPosAdcTimeWindowMax(0),
53  fHodoNegAdcTimeWindowMin(0), fHodoNegAdcTimeWindowMax(0),
54  fHodoPosInvAdcLinear(0), fHodoNegInvAdcLinear(0),
55  fHodoPosInvAdcAdc(0), fHodoNegInvAdcAdc(0), fHodoVelFit(0),
56  fHodoCableFit(0), fHodo_LCoeff(0), fHodoPos_c1(0), fHodoNeg_c1(0),
57  fHodoPos_c2(0), fHodoNeg_c2(0), fHodoSigma(0), fPosPedSum(0),
58  fPosPedSum2(0), fPosPedLimit(0), fPosPedCount(0), fNegPedSum(0),
59  fNegPedSum2(0), fNegPedLimit(0), fNegPedCount(0), fPosPed(0),
60  fPosSig(0), fPosThresh(0), fNegPed(0), fNegSig(0), fNegThresh(0)
61 {
62  // Normal constructor with name and description
63  fHodoHits = new TClonesArray("THcHodoHit",16);
64 
65  fCluster = new TClonesArray("THcScintPlaneCluster", 10);
66 
67  frPosAdcErrorFlag = new TClonesArray("THcSignalHit", 16);
68  frNegAdcErrorFlag = new TClonesArray("THcSignalHit", 16);
69 
70  frPosTDCHits = new TClonesArray("THcSignalHit",16);
71  frNegTDCHits = new TClonesArray("THcSignalHit",16);
72  frPosADCHits = new TClonesArray("THcSignalHit",16);
73  frNegADCHits = new TClonesArray("THcSignalHit",16);
74  frPosADCSums = new TClonesArray("THcSignalHit",16);
75  frNegADCSums = new TClonesArray("THcSignalHit",16);
76  frPosADCPeds = new TClonesArray("THcSignalHit",16);
77  frNegADCPeds = new TClonesArray("THcSignalHit",16);
78 
79  frPosTdcTimeRaw = new TClonesArray("THcSignalHit", 16);
80  frPosAdcPedRaw = new TClonesArray("THcSignalHit", 16);
81  frPosAdcPulseIntRaw = new TClonesArray("THcSignalHit", 16);
82  frPosAdcPulseAmpRaw = new TClonesArray("THcSignalHit", 16);
83  frPosAdcPulseTimeRaw = new TClonesArray("THcSignalHit", 16);
84 
85  frPosTdcTime = new TClonesArray("THcSignalHit", 16);
86  frPosAdcPed = new TClonesArray("THcSignalHit", 16);
87  frPosAdcPulseInt = new TClonesArray("THcSignalHit", 16);
88  frPosAdcPulseAmp = new TClonesArray("THcSignalHit", 16);
89  frPosAdcPulseTime = new TClonesArray("THcSignalHit", 16);
90 
91  frNegTdcTimeRaw = new TClonesArray("THcSignalHit", 16);
92  frNegAdcPedRaw = new TClonesArray("THcSignalHit", 16);
93  frNegAdcPulseIntRaw = new TClonesArray("THcSignalHit", 16);
94  frNegAdcPulseAmpRaw = new TClonesArray("THcSignalHit", 16);
95  frNegAdcPulseTimeRaw = new TClonesArray("THcSignalHit", 16);
96 
97  frNegTdcTime = new TClonesArray("THcSignalHit", 16);
98  frNegAdcPed = new TClonesArray("THcSignalHit", 16);
99  frNegAdcPulseInt = new TClonesArray("THcSignalHit", 16);
100  frNegAdcPulseAmp = new TClonesArray("THcSignalHit", 16);
101  frNegAdcPulseTime = new TClonesArray("THcSignalHit", 16);
102 
103 
104  fPlaneNum = planenum;
105  fTotPlanes = planenum;
106  fNScinHits = 0;
107 }
108 
109 //______________________________________________________________________________
111 {
112  // Destructor
113  if( fIsSetup )
114  RemoveVariables();
115  delete frPosAdcErrorFlag; frPosAdcErrorFlag = NULL;
116  delete frNegAdcErrorFlag; frNegAdcErrorFlag = NULL;
117 
118  delete fCluster; fCluster = NULL;
119 
120  delete fHodoHits;
121  delete frPosTDCHits;
122  delete frNegTDCHits;
123  delete frPosADCHits;
124  delete frNegADCHits;
125  delete frPosADCSums;
126  delete frNegADCSums;
127  delete frPosADCPeds;
128  delete frNegADCPeds;
129 
130  delete frPosTdcTimeRaw;
131  delete frPosAdcPedRaw;
132  delete frPosAdcPulseIntRaw;
133  delete frPosAdcPulseAmpRaw;
134  delete frPosAdcPulseTimeRaw;
135 
136  delete frPosTdcTime;
137  delete frPosAdcPed;
138  delete frPosAdcPulseInt;
139  delete frPosAdcPulseAmp;
140  delete frPosAdcPulseTime;
141 
142  delete frNegTdcTimeRaw;
143  delete frNegAdcPedRaw;
144  delete frNegAdcPulseIntRaw;
145  delete frNegAdcPulseAmpRaw;
146  delete frNegAdcPulseTimeRaw;
147 
148  delete frNegTdcTime;
149  delete frNegAdcPed;
150  delete frNegAdcPulseInt;
151  delete frNegAdcPulseAmp;
152  delete frNegAdcPulseTime;
153 
154  delete [] fPosCenter; fPosCenter = 0;
155 
156  delete [] fHodoPosMinPh; fHodoPosMinPh = NULL;
157  delete [] fHodoNegMinPh; fHodoNegMinPh = NULL;
158  delete [] fHodoPosPhcCoeff; fHodoPosPhcCoeff = NULL;
159  delete [] fHodoNegPhcCoeff; fHodoNegPhcCoeff = NULL;
160  delete [] fHodoPosTimeOffset; fHodoPosTimeOffset = NULL;
161  delete [] fHodoNegTimeOffset; fHodoNegTimeOffset = NULL;
162  delete [] fHodoPosInvAdcOffset; fHodoPosInvAdcOffset = NULL;
163  delete [] fHodoNegInvAdcOffset; fHodoNegInvAdcOffset = NULL;
164  delete [] fHodoPosInvAdcLinear; fHodoPosInvAdcLinear = NULL;
165  delete [] fHodoNegInvAdcLinear; fHodoNegInvAdcLinear = NULL;
166  delete [] fHodoPosAdcTimeWindowMax; fHodoPosAdcTimeWindowMax = NULL;
167  delete [] fHodoPosAdcTimeWindowMin; fHodoPosAdcTimeWindowMin = NULL;
168  delete [] fHodoNegAdcTimeWindowMax; fHodoNegAdcTimeWindowMax = NULL;
169  delete [] fHodoNegAdcTimeWindowMin; fHodoNegAdcTimeWindowMin = NULL;
170  delete [] fHodoNegInvAdcLinear; fHodoNegInvAdcLinear = NULL;
171  delete [] fHodoPosInvAdcAdc; fHodoPosInvAdcAdc = NULL;
172  delete [] fHodoNegInvAdcAdc; fHodoNegInvAdcAdc = NULL;
173  delete [] fHodoVelFit; fHodoVelFit = NULL;
174  delete [] fHodoCableFit; fHodoCableFit = NULL;
175  delete [] fHodo_LCoeff; fHodo_LCoeff = NULL;
176  delete [] fHodoPos_c1; fHodoPos_c1 = NULL;
177  delete [] fHodoNeg_c1; fHodoNeg_c1 = NULL;
178  delete [] fHodoPos_c2; fHodoPos_c2 = NULL;
179  delete [] fHodoNeg_c2; fHodoNeg_c2 = NULL;
180 
181 
182  delete [] fHodoVelLight; fHodoVelLight = NULL;
183  delete [] fHodoSigma; fHodoSigma = NULL;
184 
185  delete [] fPosPedSum; fPosPedSum = 0;
186  delete [] fPosPedSum2; fPosPedSum2 = 0;
187  delete [] fPosPedLimit; fPosPedLimit = 0;
188  delete [] fPosPedCount; fPosPedCount = 0;
189  delete [] fNegPedSum; fNegPedSum = 0;
190  delete [] fNegPedSum2; fNegPedSum2 = 0;
191  delete [] fNegPedLimit; fNegPedLimit = 0;
192  delete [] fNegPedCount; fNegPedCount = 0;
193  delete [] fPosPed; fPosPed = 0;
194  delete [] fNegPed; fNegPed = 0;
195  delete [] fPosThresh; fPosThresh = 0;
196  delete [] fNegThresh; fNegThresh = 0;
197 }
198 
199 //______________________________________________________________________________
200 THaAnalysisObject::EStatus THcScintillatorPlane::Init( const TDatime& date )
201 {
202 
203  // cout << "THcScintillatorPlane::Init called " << GetName() << endl;
204 
205  if( IsZombie())
206  return fStatus = kInitError;
207 
208  // How to get information for parent
209  // if( GetParent() )
210  // fOrigin = GetParent()->GetOrigin();
211 
212  EStatus status;
213  if( (status=THaSubDetector::Init( date )) )
214  return fStatus = status;
215 
216  // Get the Hodoscope hitlist
217  // Can't seem to cast to THcHitList. What to do if we want to use
218  // THcScintillatorPlane as a subdetector to other than THcHodoscope?
219  // fParentHitList = static_cast<THcHodoscope*>(GetParent())->GetHitList();
220 
221  return fStatus = kOK;
222 
223 }
224 
225 //_____________________________________________________________________________
227 {
228 
229  // See what file it looks for
230 
231  // static const char* const here = "ReadDatabase()";
232  char prefix[2];
233 
234  prefix[0]=tolower(GetParent()->GetPrefix()[0]);
235  prefix[1]='\0';
236 
237  // need this further down so read them first! GN
238  string parname = "scin_" + string(GetName()) + "_nr";
239  DBRequest list_1[] = {
240  {parname.c_str(), &fNelem, kInt},
241  {0}
242  };
243  gHcParms->LoadParmValues(list_1, prefix);
244 
245  // Based on the signs of these quantities in the .pos file the correspondence
246  // should be bot=>left and top=>right when comparing x and y-type scintillators
247  char tmpleft[6], tmpright[6];
248  if (fPlaneNum==1 || fPlaneNum==3) {
249  strcpy(tmpleft,"left");
250  strcpy(tmpright,"right");
251  }
252  else {
253  strcpy(tmpleft,"bot");
254  strcpy(tmpright,"top");
255  }
256 
257  delete [] fPosCenter; fPosCenter = new Double_t[fNelem];
258 
259  DBRequest list[]={
260  {Form("scin_%s_zpos",GetName()), &fZpos, kDouble},
261  {Form("scin_%s_dzpos",GetName()), &fDzpos, kDouble},
262  {Form("scin_%s_size",GetName()), &fSize, kDouble},
263  {Form("scin_%s_spacing",GetName()), &fSpacing, kDouble},
264  {Form("scin_%s_%s",GetName(),tmpleft), &fPosLeft,kDouble},
265  {Form("scin_%s_%s",GetName(),tmpright), &fPosRight,kDouble},
266  {Form("scin_%s_offset",GetName()), &fPosOffset, kDouble},
267  {Form("scin_%s_center",GetName()), fPosCenter,kDouble,fNelem},
268  {"tofusinginvadc", &fTofUsingInvAdc, kInt, 0, 1},
269  {"hodo_adc_mode", &fADCMode, kInt, 0, 1},
270  {"hodo_pedestal_scale", &fADCPedScaleFactor, kDouble, 0, 1},
271  {"hodo_adc_diag_cut", &fADCDiagCut, kInt, 0, 1},
272  {"cosmicflag", &fCosmicFlag, kInt, 0, 1},
273  {"hodo_debug_adc", &fDebugAdc, kInt, 0, 1},
274  {0}
275  };
276 
277  fDebugAdc = 0; // Set ADC debug parameter to false unless set in parameter file
278  fTofUsingInvAdc = 1;
279  //fADCMode = kADCStandard;
280  fADCMode = kADCDynamicPedestal;
281  fADCPedScaleFactor = 1.0;
282  fADCDiagCut = 50.0;
283  fCosmicFlag=0;
284  gHcParms->LoadParmValues((DBRequest*)&list,prefix);
285  if (fCosmicFlag==1) cout << " setup for cosmics in scint plane"<< endl;
286  // cout << " cosmic flag = " << fCosmicFlag << endl;
287  // fetch the parameter from the temporary list
288 
289  // Retrieve parameters we need from parent class
290  // Common for all planes
291 
292  THcHodoscope* parent = (THcHodoscope*)GetParent();
293  fHodoSlop= parent->GetHodoSlop(fPlaneNum-1);
294  fTdcOffset= parent->GetTdcOffset(fPlaneNum-1);
295  fAdcTdcOffset= parent->GetAdcTdcOffset(fPlaneNum-1);
296  fScinTdcMin=parent->GetTdcMin();
297  fScinTdcMax=parent->GetTdcMax();
298  fScinTdcToTime=parent->GetTdcToTime();
299  fTofTolerance=parent->GetTofTolerance();
300  fBetaNominal=parent->GetBetaNominal();
301  fStartTimeCenter=parent->GetStartTimeCenter();
302  fStartTimeSlop=parent->GetStartTimeSlop();
303  // Parameters for this plane
304  fHodoPosMinPh = new Double_t[fNelem];
305  fHodoNegMinPh = new Double_t[fNelem];
306  fHodoPosPhcCoeff = new Double_t[fNelem];
307  fHodoNegPhcCoeff = new Double_t[fNelem];
308  fHodoPosTimeOffset = new Double_t[fNelem];
309  fHodoNegTimeOffset = new Double_t[fNelem];
310  fHodoVelLight = new Double_t[fNelem];
311  fHodoPosInvAdcOffset = new Double_t[fNelem];
312  fHodoNegInvAdcOffset = new Double_t[fNelem];
313  fHodoPosInvAdcLinear = new Double_t[fNelem];
314  fHodoNegInvAdcLinear = new Double_t[fNelem];
315  fHodoPosAdcTimeWindowMin = new Double_t[fNelem];
316  fHodoNegAdcTimeWindowMin = new Double_t[fNelem];
317  fHodoPosAdcTimeWindowMax = new Double_t[fNelem];
318  fHodoNegAdcTimeWindowMax = new Double_t[fNelem];
319  fHodoPosInvAdcAdc = new Double_t[fNelem];
320  fHodoNegInvAdcAdc = new Double_t[fNelem];
321  fHodoSigma = new Double_t[fNelem];
322 
323  //New Time-Walk Calibration Parameters
324  fHodoVelFit=new Double_t [fNelem];
325  fHodoCableFit=new Double_t [fNelem];
326  fHodo_LCoeff=new Double_t [fNelem];
327  fHodoPos_c1=new Double_t [fNelem];
328  fHodoNeg_c1=new Double_t [fNelem];
329  fHodoPos_c2=new Double_t [fNelem];
330  fHodoNeg_c2=new Double_t [fNelem];
331 
332  for(Int_t j=0;j<(Int_t) fNelem;j++) {
333  Int_t index=parent->GetScinIndex(fPlaneNum-1,j);
334  fHodoPosAdcTimeWindowMin[j] = parent->GetHodoPosAdcTimeWindowMin(index);
335  fHodoPosAdcTimeWindowMax[j] = parent->GetHodoPosAdcTimeWindowMax(index);
336  fHodoNegAdcTimeWindowMin[j] = parent->GetHodoNegAdcTimeWindowMin(index);
337  fHodoNegAdcTimeWindowMax[j] = parent->GetHodoNegAdcTimeWindowMax(index);
338  fHodoPosMinPh[j] = parent->GetHodoPosMinPh(index);
339  fHodoNegMinPh[j] = parent->GetHodoNegMinPh(index);
340  fHodoPosPhcCoeff[j] = parent->GetHodoPosPhcCoeff(index);
341  fHodoNegPhcCoeff[j] = parent->GetHodoNegPhcCoeff(index);
342  fHodoPosTimeOffset[j] = parent->GetHodoPosTimeOffset(index);
343  fHodoNegTimeOffset[j] = parent->GetHodoNegTimeOffset(index);
344  fHodoPosInvAdcOffset[j] = parent->GetHodoPosInvAdcOffset(index);
345  fHodoNegInvAdcOffset[j] = parent->GetHodoNegInvAdcOffset(index);
346  fHodoPosInvAdcLinear[j] = parent->GetHodoPosInvAdcLinear(index);
347  fHodoNegInvAdcLinear[j] = parent->GetHodoNegInvAdcLinear(index);
348  fHodoPosInvAdcAdc[j] = parent->GetHodoPosInvAdcAdc(index);
349  fHodoNegInvAdcAdc[j] = parent->GetHodoNegInvAdcAdc(index);
350  fHodoVelLight[j] = parent->GetHodoVelLight(index);
351  //Get Time-Walk correction param
352  fHodoVelFit[j] = parent->GetHodoVelFit(index);
353  fHodoCableFit[j] = parent->GetHodoCableFit(index);
354  fHodo_LCoeff[j] = parent->GetHodoLCoeff(index);
355  fHodoPos_c1[j] = parent->GetHodoPos_c1(index);
356  fHodoNeg_c1[j] = parent->GetHodoNeg_c1(index);
357  fHodoPos_c2[j] = parent->GetHodoPos_c2(index);
358  fHodoNeg_c2[j] = parent->GetHodoNeg_c2(index);
359 
360  Double_t possigma = parent->GetHodoPosSigma(index);
361  Double_t negsigma = parent->GetHodoNegSigma(index);
362  fHodoSigma[j] = TMath::Sqrt(possigma*possigma+negsigma*negsigma)/2.0;
363 
364 
365 
366 
367  }
368 
369  fTdc_Thrs = parent->GetTDCThrs();
370  // cout <<" plane num = "<<fPlaneNum<<endl;
371  // cout <<" nelem = "<<fNelem<<endl;
372  // cout <<" zpos = "<<fZpos<<endl;
373  // cout <<" dzpos = "<<fDzpos<<endl;
374  // cout <<" spacing = "<<fSpacing<<endl;
375  // cout <<" size = "<<fSize<<endl;
376  // cout <<" hodoslop = "<<fHodoSlop<<endl;
377  // cout <<"PosLeft = "<<fPosLeft<<endl;
378  // cout <<"PosRight = "<<fPosRight<<endl;
379  // cout <<"PosOffset = "<<fPosOffset<<endl;
380  // cout <<"PosCenter[0] = "<<fPosCenter[0]<<endl;
381 
382  // Think we will make special methods to pass most
383  // How generic do we want to make this class?
384  // The way we get parameter data is going to be pretty specific to
385  // our parameter file naming conventions. But on the other hand,
386  // the Hall A analyzer read database is pretty specific.
387  // Is there any way for this class to know which spectrometer he
388  // belongs too?
389 
390  // Create arrays to hold results here
391  InitializePedestals();
392 
393 
394  fNumGoodPosAdcHits = vector<Int_t> (fNelem, 0.0);
395  fNumGoodNegAdcHits = vector<Int_t> (fNelem, 0.0);
396  fNumGoodPosTdcHits = vector<Int_t> (fNelem, 0.0);
397  fNumGoodNegTdcHits = vector<Int_t> (fNelem, 0.0);
398 
399  fGoodPosAdcPed = vector<Double_t> (fNelem, 0.0);
400  fGoodNegAdcPed = vector<Double_t> (fNelem, 0.0);
401  fGoodPosAdcMult = vector<Double_t> (fNelem, 0.0);
402  fGoodNegAdcMult = vector<Double_t> (fNelem, 0.0);
403  fGoodPosAdcHitUsed = vector<Double_t> (fNelem, 0.0);
404  fGoodNegAdcHitUsed = vector<Double_t> (fNelem, 0.0);
405  fGoodPosAdcPulseAmp = vector<Double_t> (fNelem, 0.0);
406  fGoodNegAdcPulseAmp = vector<Double_t> (fNelem, 0.0);
407  fGoodPosAdcPulseInt = vector<Double_t> (fNelem, 0.0);
408  fGoodNegAdcPulseInt = vector<Double_t> (fNelem, 0.0);
409  fGoodPosAdcPulseTime = vector<Double_t> (fNelem, 0.0);
410  fGoodNegAdcPulseTime = vector<Double_t> (fNelem, 0.0);
411  fGoodPosAdcTdcDiffTime = vector<Double_t> (fNelem, 0.0);
412  fGoodNegAdcTdcDiffTime = vector<Double_t> (fNelem, 0.0);
413 
414  fGoodPosTdcTimeUnCorr = vector<Double_t> (fNelem, 0.0);
415  fGoodNegTdcTimeUnCorr = vector<Double_t> (fNelem, 0.0);
416  fGoodPosTdcTimeCorr = vector<Double_t> (fNelem, 0.0);
417  fGoodNegTdcTimeCorr = vector<Double_t> (fNelem, 0.0);
418  fGoodPosTdcTimeTOFCorr = vector<Double_t> (fNelem, 0.0);
419  fGoodNegTdcTimeTOFCorr = vector<Double_t> (fNelem, 0.0);
420  fGoodPosTdcTimeWalkCorr = vector<Double_t> (fNelem, 0.0);
421  fGoodNegTdcTimeWalkCorr = vector<Double_t> (fNelem, 0.0);
422  fGoodDiffDistTrack = vector<Double_t> (fNelem, 0.0);
423  return kOK;
424 }
425 //_____________________________________________________________________________
427 {
432  // cout << "THcScintillatorPlane::DefineVariables called " << GetName() << endl;
433 
434  if( mode == kDefine && fIsSetup ) return kOK;
435  fIsSetup = ( mode == kDefine );
436 
437  // Register variables in global list
438 
439  if (fDebugAdc) {
440  RVarDef vars[] = {
441  {"posAdcErrorFlag", "Error Flag for When FPGA Fails", "frPosAdcErrorFlag.THcSignalHit.GetData()"},
442  {"negAdcErrorFlag", "Error Flag for When FPGA Fails", "frNegAdcErrorFlag.THcSignalHit.GetData()"},
443 
444  {"posTdcTimeRaw", "List of positive raw TDC values.", "frPosTdcTimeRaw.THcSignalHit.GetData()"},
445  {"posAdcPedRaw", "List of positive raw ADC pedestals", "frPosAdcPedRaw.THcSignalHit.GetData()"},
446  {"posAdcPulseIntRaw", "List of positive raw ADC pulse integrals.", "frPosAdcPulseIntRaw.THcSignalHit.GetData()"},
447  {"posAdcPulseAmpRaw", "List of positive raw ADC pulse amplitudes.", "frPosAdcPulseAmpRaw.THcSignalHit.GetData()"},
448  {"posAdcPulseTimeRaw", "List of positive raw ADC pulse times.", "frPosAdcPulseTimeRaw.THcSignalHit.GetData()"},
449 
450  {"posTdcTime", "List of positive TDC values.", "frPosTdcTime.THcSignalHit.GetData()"},
451  {"posAdcPed", "List of positive ADC pedestals", "frPosAdcPed.THcSignalHit.GetData()"},
452  {"posAdcPulseInt", "List of positive ADC pulse integrals.", "frPosAdcPulseInt.THcSignalHit.GetData()"},
453  {"posAdcPulseAmp", "List of positive ADC pulse amplitudes.", "frPosAdcPulseAmp.THcSignalHit.GetData()"},
454  {"posAdcPulseTime", "List of positive ADC pulse times.", "frPosAdcPulseTime.THcSignalHit.GetData()"},
455 
456  {"negTdcTimeRaw", "List of negative raw TDC values.", "frNegTdcTimeRaw.THcSignalHit.GetData()"},
457  {"negAdcPedRaw", "List of negative raw ADC pedestals", "frNegAdcPedRaw.THcSignalHit.GetData()"},
458  {"negAdcPulseIntRaw", "List of negative raw ADC pulse integrals.", "frNegAdcPulseIntRaw.THcSignalHit.GetData()"},
459  {"negAdcPulseAmpRaw", "List of negative raw ADC pulse amplitudes.", "frNegAdcPulseAmpRaw.THcSignalHit.GetData()"},
460  {"negAdcPulseTimeRaw", "List of negative raw ADC pulse times.", "frNegAdcPulseTimeRaw.THcSignalHit.GetData()"},
461 
462  {"negTdcTime", "List of negative TDC values.", "frNegTdcTime.THcSignalHit.GetData()"},
463  {"negAdcPed", "List of negative ADC pedestals", "frNegAdcPed.THcSignalHit.GetData()"},
464  {"negAdcPulseInt", "List of negative ADC pulse integrals.", "frNegAdcPulseInt.THcSignalHit.GetData()"},
465  {"negAdcPulseAmp", "List of negative ADC pulse amplitudes.", "frNegAdcPulseAmp.THcSignalHit.GetData()"},
466  {"negAdcPulseTime", "List of negative ADC pulse times.", "frNegAdcPulseTime.THcSignalHit.GetData()"},
467 
468  {"totNumPosAdcHits", "Total Number of Positive ADC Hits", "fTotNumPosAdcHits"}, // Hodo+ raw ADC multiplicity Int_t
469  {"totNumNegAdcHits", "Total Number of Negative ADC Hits", "fTotNumNegAdcHits"}, // Hodo- raw ADC multiplicity ""
470  {"totNumAdcHits", "Total Number of PMTs Hit (as measured by ADCs)", "fTotNumAdcHits"}, // Hodo raw ADC multiplicity ""
471 
472  {"totNumPosTdcHits", "Total Number of Positive TDC Hits", "fTotNumPosTdcHits"}, // Hodo+ raw TDC multiplicity ""
473  {"totNumNegTdcHits", "Total Number of Negative TDC Hits", "fTotNumNegTdcHits"}, // Hodo- raw TDC multiplicity ""
474  {"totNumTdcHits", "Total Number of PMTs Hits (as measured by TDCs)", "fTotNumTdcHits"}, // Hodo raw TDC multiplicity ""
475  { 0 }
476  };
477  DefineVarsFromList( vars, mode);
478  } //end debug statement
479 
480  RVarDef vars[] = {
481  {"nhits", "Number of paddle hits (passed TDC && ADC Min and Max cuts for either end)", "GetNScinHits() "},
482 
483  {"posTdcCounter", "List of positive TDC counter numbers.", "frPosTdcTimeRaw.THcSignalHit.GetPaddleNumber()"}, //Hodo+ raw TDC occupancy
484  {"posAdcCounter", "List of positive ADC counter numbers.", "frPosAdcPulseIntRaw.THcSignalHit.GetPaddleNumber()"}, //Hodo+ raw ADC occupancy
485  {"negTdcCounter", "List of negative TDC counter numbers.", "frNegTdcTimeRaw.THcSignalHit.GetPaddleNumber()"}, //Hodo- raw TDC occupancy
486  {"negAdcCounter", "List of negative ADC counter numbers.", "frNegAdcPulseIntRaw.THcSignalHit.GetPaddleNumber()"}, //Hodo- raw ADC occupancy
487 
488  {"fptime", "Time at focal plane", "GetFpTime()"},
489 
490  {"numGoodPosAdcHits", "Number of Good Positive ADC Hits Per PMT", "fNumGoodPosAdcHits"}, // Hodo+ good ADC occupancy - vector<Int_t>
491  {"numGoodNegAdcHits", "Number of Good Negative ADC Hits Per PMT", "fNumGoodNegAdcHits"}, // Hodo- good ADC occupancy - vector <Int_t>
492 
493  {"numGoodPosTdcHits", "Number of Good Positive TDC Hits Per PMT", "fNumGoodPosTdcHits"}, // Hodo+ good TDC occupancy - vector<Int_t>
494  {"numGoodNegTdcHits", "Number of Good Negative TDC Hits Per PMT", "fNumGoodNegTdcHits"}, // Hodo- good TDC occupancy - vector <Int_t>
495 
496 
497  {"totNumGoodPosAdcHits", "Total Number of Good Positive ADC Hits", "fTotNumGoodPosAdcHits"}, // Hodo+ good ADC multiplicity - Int_t
498  {"totNumGoodNegAdcHits", "Total Number of Good Negative ADC Hits", "fTotNumGoodNegAdcHits"}, // Hodo- good ADC multiplicity - Int_t
499  {"totNumGoodAdcHits", "TotalNumber of Good ADC Hits Per PMT", "fTotNumGoodAdcHits"}, // Hodo good ADC multiplicity - Int_t
500 
501  {"totNumGoodPosTdcHits", "Total Number of Good Positive TDC Hits", "fTotNumGoodPosTdcHits"}, // Hodo+ good TDC multiplicity - Int_t
502  {"totNumGoodNegTdcHits", "Total Number of Good Negative TDC Hits", "fTotNumGoodNegTdcHits"}, // Hodo- good TDC multiplicity - Int_t
503  {"totNumGoodTdcHits", "TotalNumber of Good TDC Hits Per PMT", "fTotNumGoodTdcHits"}, // Hodo good TDC multiplicity - Int_t
504 
505 
506 
507  // {"GoodPaddle", "List of Paddle Numbers (passed TDC && ADC Min and Max cuts for either end)", "fHodoHits.THcHodoHit.GetPaddleNumber()"},
508 
509  {"GoodPosAdcPed", "List of Positive ADC pedestals (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcPed"}, //vector<Double_t>
510  {"GoodNegAdcPed", "List of Negative ADC pedestals (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcPed"}, //vector<Double_t>
511  {"GoodPosAdcMult", "List of Positive ADC Mult (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcMult"}, //vector<Double_t>
512  {"GoodNegAdcMult", "List of Negative ADC Mult (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcMult"}, //vector<Double_t>
513  {"GoodPosAdcHitUsed", "List of Positive ADC Hit Used (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcHitUsed"}, //vector<Double_t>
514  {"GoodNegAdcHitUsed", "List of Negative ADC Hit Used (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcHitUsed"}, //vector<Double_t>
515 
516  {"GoodNegTdcTimeUnCorr", "List of negative TDC values (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegTdcTimeUnCorr"}, //Units ns
517  {"GoodNegTdcTimeCorr", "List of negative corrected TDC values (corrected for PMT offset and ADC)", "fGoodNegTdcTimeCorr"},
518  {"GoodNegTdcTimeTOFCorr", "List of negative corrected TDC values (corrected for TOF)", "fGoodNegTdcTimeTOFCorr"},
519  {"GoodNegTdcTimeWalkCorr", "List of negative corrected TDC values (corrected for Time-Walk)", "fGoodNegTdcTimeWalkCorr"},
520  {"GoodNegAdcPulseInt", "List of negative ADC values (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcPulseInt"},
521  {"GoodPosTdcTimeUnCorr", "List of positive TDC values (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosTdcTimeUnCorr"},
522  {"GoodPosTdcTimeCorr", "List of positive corrected TDC values (corrected for PMT offset and ADC)", "fGoodPosTdcTimeCorr"},
523  {"GoodPosTdcTimeTOFCorr", "List of positive corrected TDC values (corrected for TOF)", "fGoodPosTdcTimeTOFCorr"},
524  {"GoodPosTdcTimeWalkCorr", "List of positive corrected TDC values (corrected for Time-Walk)", "fGoodPosTdcTimeWalkCorr"},
525  {"GoodPosAdcPulseInt", "List of positive ADC values (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcPulseInt"},
526  {"GoodPosAdcPulseAmp", "List of positive ADC peak amp (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcPulseAmp"},
527  {"GoodNegAdcPulseAmp", "List of Negative ADC peak amp (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcPulseAmp"},
528  {"GoodPosAdcPulseTime", "List of positive ADC time (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcPulseTime"},
529  {"GoodNegAdcPulseTime", "List of Negative ADC time (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcPulseTime"},
530  {"GoodPosAdcTdcDiffTime", "List of positive TDC - ADC time (passed TDC && ADC Min and Max cuts for either end)", "fGoodPosAdcTdcDiffTime"},
531  {"GoodNegAdcTdcDiffTime", "List of Negative TDC - ADC time (passed TDC && ADC Min and Max cuts for either end)", "fGoodNegAdcTdcDiffTime"},
532  {"DiffDisTrack", "Difference between track and scintillator position (cm)", "fHitDistance"},
533  {"DiffDisTrackCorr", "TW Corrected Dist Difference between track and scintillator position (cm)", "fGoodDiffDistTrack"},
534  {"TrackXPos", "Track X position at plane (cm)", "fTrackXPosition"},
535  {"TrackYPos", "Track Y position at plane (cm)", "fTrackYPosition"},
536  {"ScinXPos", "Scint Average Y position at plane (cm)", "fScinXPos"},
537  {"ScinYPos", "Scint Average Xposition at plane (cm)", "fScinYPos"},
538  {"NumClus", "Number of clusters", "fNumberClusters"},
539  {"Clus.Pos", "Position of each paddle clusters", "fCluster.THcScintPlaneCluster.GetClusterPosition()"},
540  {"Clus.Size", "Size of each paddle clusters", "fCluster.THcScintPlaneCluster.GetClusterSize()"},
541  {"Clus.Flag", "Flag of each paddle clusters", "fCluster.THcScintPlaneCluster.GetClusterFlag()"},
542  {"Clus.UsedFlag", "USed Flag of each paddle clusters", "fCluster.THcScintPlaneCluster.GetClusterUsedFlag()"},
543  {"PosTdcRefTime", "Reference time of Pos TDC", "fPosTdcRefTime"},
544  {"NegTdcRefTime", "Reference time of Neg TDC", "fNegTdcRefTime"},
545  {"PosAdcRefTime", "Reference time of Pos ADC", "fPosAdcRefTime"},
546  {"NegAdcRefTime", "Reference time of Neg aDC", "fNegAdcRefTime"},
547  {"PosTdcRefDiffTime", "Reference Diff time of Pos TDC", "fPosTdcRefDiffTime"},
548  {"NegTdcRefDiffTime", "Reference Diff time of Neg TDC", "fNegTdcRefDiffTime"},
549  {"PosAdcRefDiffTime", "Reference Diff time of Pos ADC", "fPosAdcRefDiffTime"},
550  {"NegAdcRefDiffTime", "Reference Diff time of Neg aDC", "fNegAdcRefDiffTime"},
551  //{"ngoodhits", "Number of paddle hits (passed tof tolerance and used to determine the focal plane time )", "GetNGoodHits() "},
552  { 0 }
553  };
554 
555  return DefineVarsFromList(vars, mode);
556 
557 }
558 //_____________________________________________________________________________
560 {
565  //cout << " Calling THcScintillatorPlane::Clear " << GetName() << endl;
566  // Clears the hit lists
567  fCluster->Clear();
568 
569  frPosAdcErrorFlag->Clear();
570  frNegAdcErrorFlag->Clear();
571 
572  fHodoHits->Clear();
573  frPosTDCHits->Clear();
574  frNegTDCHits->Clear();
575  frPosADCHits->Clear();
576  frNegADCHits->Clear();
577 
578  frPosTdcTimeRaw->Clear();
579  frPosAdcPedRaw->Clear();
580  frPosAdcPulseIntRaw->Clear();
581  frPosAdcPulseAmpRaw->Clear();
582  frPosAdcPulseTimeRaw->Clear();
583 
584  frPosTdcTime->Clear();
585  frPosAdcPed->Clear();
586  frPosAdcPulseInt->Clear();
587  frPosAdcPulseAmp->Clear();
588  frPosAdcPulseTime->Clear();
589 
590  frNegTdcTimeRaw->Clear();
591  frNegAdcPedRaw->Clear();
592  frNegAdcPulseIntRaw->Clear();
593  frNegAdcPulseAmpRaw->Clear();
594  frNegAdcPulseTimeRaw->Clear();
595 
596  frNegTdcTime->Clear();
597  frNegAdcPed->Clear();
598  frNegAdcPulseInt->Clear();
599  frNegAdcPulseAmp->Clear();
600  frNegAdcPulseTime->Clear();
601 
602 
603  //Clear occupancies
604  for (UInt_t ielem = 0; ielem < fNumGoodPosAdcHits.size(); ielem++)
605  fNumGoodPosAdcHits.at(ielem) = 0;
606  for (UInt_t ielem = 0; ielem < fNumGoodNegAdcHits.size(); ielem++)
607  fNumGoodNegAdcHits.at(ielem) = 0;
608 
609  for (UInt_t ielem = 0; ielem < fNumGoodPosTdcHits.size(); ielem++)
610  fNumGoodPosTdcHits.at(ielem) = 0;
611  for (UInt_t ielem = 0; ielem < fNumGoodNegTdcHits.size(); ielem++)
612  fNumGoodNegTdcHits.at(ielem) = 0;
613 
614  //Clear Ped/Amps/Int/Time
615  for (UInt_t ielem = 0; ielem < fGoodPosAdcPed.size(); ielem++) {
616  fGoodPosAdcPed.at(ielem) = 0.0;
617  fGoodPosAdcMult.at(ielem) = 0.0;
618  fGoodPosAdcHitUsed.at(ielem) = 0.0;
619  fGoodPosAdcPulseInt.at(ielem) = 0.0;
620  fGoodPosAdcPulseAmp.at(ielem) = 0.0;
621  fGoodPosAdcPulseTime.at(ielem) = kBig;
622  fGoodPosAdcTdcDiffTime.at(ielem) = kBig;
623  }
624  for (UInt_t ielem = 0; ielem < fGoodNegAdcPed.size(); ielem++) {
625  fGoodNegAdcPed.at(ielem) = 0.0;
626  fGoodNegAdcMult.at(ielem) = 0.0;
627  fGoodNegAdcHitUsed.at(ielem) = 0.0;
628  fGoodNegAdcPulseInt.at(ielem) = 0.0;
629  fGoodNegAdcPulseAmp.at(ielem) = 0.0;
630  fGoodNegAdcPulseTime.at(ielem) = kBig;
631  fGoodNegAdcTdcDiffTime.at(ielem) = kBig;
632  }
633 
634  //Clear Good TDC Variables
635  for (UInt_t ielem = 0; ielem < fGoodPosTdcTimeUnCorr.size(); ielem++) {
636  fGoodPosTdcTimeUnCorr.at(ielem) = kBig;
637  fGoodPosTdcTimeCorr.at(ielem) = kBig;
638  fGoodPosTdcTimeTOFCorr.at(ielem) = kBig;
639  fGoodPosTdcTimeWalkCorr.at(ielem) = kBig;
640  }
641 
642  for (UInt_t ielem = 0; ielem < fGoodNegTdcTimeUnCorr.size(); ielem++) {
643  fGoodNegTdcTimeUnCorr.at(ielem) = kBig;
644  fGoodNegTdcTimeCorr.at(ielem) = kBig;
645  fGoodNegTdcTimeTOFCorr.at(ielem) = kBig;
646  fGoodNegTdcTimeWalkCorr.at(ielem) = kBig;
647  }
648 
649  for (UInt_t ielem = 0; ielem < fGoodDiffDistTrack.size(); ielem++) {
650  fGoodDiffDistTrack.at(ielem) = kBig;
651  }
652 
653  fpTime = -1.e4;
654  fHitDistance = kBig;
655  fScinYPos = kBig;
656  fScinXPos = kBig;
657  fTrackXPosition = kBig;
658  fTrackYPosition = kBig;
659  fNumberClusters=0;
660  fPosTdcRefTime = kBig;
661  fPosAdcRefTime = kBig;
662  fNegTdcRefTime = kBig;
663  fNegAdcRefTime = kBig;
664  fPosTdcRefDiffTime = kBig;
665  fPosAdcRefDiffTime = kBig;
666  fNegTdcRefDiffTime = kBig;
667  fNegAdcRefDiffTime = kBig;
668 }
669 
670 //_____________________________________________________________________________
672 {
675  cout << " Calling THcScintillatorPlane::Decode " << GetName() << endl;
676 
677  return 0;
678 }
679 //_____________________________________________________________________________
681 {
685  cout <<"*******************************\n";
686  cout <<"NOW IN THcScintilatorPlane::CoarseProcess!!!!\n";
687  cout <<"*******************************\n";
688  // HitCount();
689 
690  return 0;
691 }
692 
693 //_____________________________________________________________________________
695 {
698  return 0;
699 }
700 
701 
702 //_____________________________________________________________________________
704 {
721  //raw
722  fPosTdcRefTime = kBig;
723  fPosAdcRefTime = kBig;
724  fNegTdcRefTime = kBig;
725  fNegAdcRefTime = kBig;
726  fPosTdcRefDiffTime = kBig;
727  fPosAdcRefDiffTime = kBig;
728  fNegTdcRefDiffTime = kBig;
729  fNegAdcRefDiffTime = kBig;
730  Int_t nrPosTDCHits=0;
731  Int_t nrNegTDCHits=0;
732  Int_t nrPosADCHits=0;
733  Int_t nrNegADCHits=0;
734  UInt_t nrPosAdcHits = 0;
735  UInt_t nrNegAdcHits = 0;
736  UInt_t nrPosTdcHits = 0;
737  UInt_t nrNegTdcHits = 0;
738  frPosTDCHits->Clear();
739  frNegTDCHits->Clear();
740  frPosADCHits->Clear();
741  frNegADCHits->Clear();
742  frPosADCSums->Clear();
743  frNegADCSums->Clear();
744  frPosADCPeds->Clear();
745  frNegADCPeds->Clear();
746 
747 
748  frPosTdcTimeRaw->Clear();
749  frPosAdcPedRaw->Clear();
750  frPosAdcPulseIntRaw->Clear();
751  frPosAdcPulseAmpRaw->Clear();
752  frPosAdcPulseTimeRaw->Clear();
753 
754  frPosTdcTime->Clear();
755  frPosAdcPed->Clear();
756  frPosAdcPulseInt->Clear();
757  frPosAdcPulseAmp->Clear();
758  frPosAdcPulseTime->Clear();
759 
760  frNegTdcTimeRaw->Clear();
761  frNegAdcPedRaw->Clear();
762  frNegAdcPulseIntRaw->Clear();
763  frNegAdcPulseAmpRaw->Clear();
764  frNegAdcPulseTimeRaw->Clear();
765 
766  frNegTdcTime->Clear();
767  frNegAdcPed->Clear();
768  frNegAdcPulseInt->Clear();
769  frNegAdcPulseAmp->Clear();
770  frNegAdcPulseTime->Clear();
771  //stripped
772  fNScinHits=0;
773 
774  fTotNumGoodPosAdcHits = 0;
775  fTotNumGoodNegAdcHits = 0;
776  fTotNumGoodAdcHits = 0;
777 
778  fTotNumPosAdcHits = 0;
779  fTotNumNegAdcHits = 0;
780  fTotNumAdcHits = 0;
781 
782 
783  fTotNumPosTdcHits = 0;
784  fTotNumNegTdcHits = 0;
785  fTotNumTdcHits = 0;
786 
787  fTotNumGoodPosTdcHits = 0;
788  fTotNumGoodNegTdcHits = 0;
789  fTotNumGoodTdcHits = 0;
790 
791 
792  fHodoHits->Clear();
793  Int_t nrawhits = rawhits->GetLast()+1;
794  // cout << "THcScintillatorPlane::ProcessHits " << fPlaneNum << " " << nexthit << "/" << nrawhits << endl;
795  Int_t ihit = nexthit;
796 
797  //cout << "THcScintillatorPlane: " << GetName() << " raw hits = " << nrawhits << endl;
798 
799  // A THcRawHodoHit contains all the information (tdc and adc for both
800  // pmts) for a single paddle for a single trigger. The tdc information
801  // might include multiple hits if it uses a multihit tdc.
802  // Use "ihit" as the index over THcRawHodoHit objects. Use
803  // "thit" to index over multiple tdc hits within an "ihit".
804  Bool_t problem_flag=kFALSE;
805  while(ihit < nrawhits) {
806  THcRawHodoHit* hit = (THcRawHodoHit *) rawhits->At(ihit);
807  if(hit->fPlane > fPlaneNum) {
808  break;
809  }
810  Int_t padnum=hit->fCounter;
811 
812  Int_t index=padnum-1;
813 
814 
815  THcRawTdcHit& rawPosTdcHit = hit->GetRawTdcHitPos();
816  if (rawPosTdcHit.GetNHits() >0 && rawPosTdcHit.HasRefTime()) {
817  if (fPosTdcRefTime == kBig) {
818  fPosTdcRefTime=rawPosTdcHit.GetRefTime() ;
819  fPosTdcRefDiffTime=rawPosTdcHit.GetRefDiffTime() ;
820  }
821  if (fPosTdcRefTime != rawPosTdcHit.GetRefTime()) {
822  cout << "THcScintillatorPlane: " << GetName() << " reftime problem at paddle num = " << padnum << " TDC pos hits = " << rawPosTdcHit.GetNHits() << endl;
823  problem_flag=kTRUE;
824  }
825  }
826  for (UInt_t thit=0; thit<rawPosTdcHit.GetNHits(); ++thit) {
827  ((THcSignalHit*) frPosTdcTimeRaw->ConstructedAt(nrPosTdcHits))->Set(padnum, rawPosTdcHit.GetTimeRaw(thit));
828  ((THcSignalHit*) frPosTdcTime->ConstructedAt(nrPosTdcHits))->Set(padnum, rawPosTdcHit.GetTime(thit));
829  ++nrPosTdcHits;
830  fTotNumTdcHits++;
831  fTotNumPosTdcHits++;
832  }
833  THcRawTdcHit& rawNegTdcHit = hit->GetRawTdcHitNeg();
834  if (rawNegTdcHit.GetNHits() >0 && rawNegTdcHit.HasRefTime()) {
835  if (fNegTdcRefTime == kBig) {
836  fNegTdcRefTime=rawNegTdcHit.GetRefTime() ;
837  fNegTdcRefDiffTime=rawNegTdcHit.GetRefDiffTime() ;
838  }
839  if (fNegTdcRefTime != rawNegTdcHit.GetRefTime()) {
840  cout << "THcScintillatorPlane: " << GetName()<< " Neg TDC reftime problem at paddle num = " << padnum << " TDC neg hits = " << rawNegTdcHit.GetNHits() << endl;
841  problem_flag=kTRUE;
842  }
843  }
844  // cout << " paddle num = " << padnum << " TDC Neg hits = " << rawNegTdcHit.GetNHits() << endl;
845  for (UInt_t thit=0; thit<rawNegTdcHit.GetNHits(); ++thit) {
846  ((THcSignalHit*) frNegTdcTimeRaw->ConstructedAt(nrNegTdcHits))->Set(padnum, rawNegTdcHit.GetTimeRaw(thit));
847  ((THcSignalHit*) frNegTdcTime->ConstructedAt(nrNegTdcHits))->Set(padnum, rawNegTdcHit.GetTime(thit));
848  ++nrNegTdcHits;
849  fTotNumTdcHits++;
850  fTotNumNegTdcHits++;
851  }
852  THcRawAdcHit& rawPosAdcHit = hit->GetRawAdcHitPos();
853  if (rawPosAdcHit.GetNPulses() >0 && rawPosAdcHit.HasRefTime()) {
854  if (fPosAdcRefTime == kBig ) {
855  fPosAdcRefTime=rawPosAdcHit.GetRefTime() ;
856  fPosAdcRefDiffTime=rawPosAdcHit.GetRefDiffTime() ;
857  }
858  if (fPosAdcRefTime != rawPosAdcHit.GetRefTime()) {
859  cout << "THcScintillatorPlane: " << GetName()<< " Pos ADC reftime problem at paddle num = " << padnum << " ADC pos hits = " << rawPosAdcHit.GetNPulses() << endl;
860  problem_flag=kTRUE;
861  }
862  }
863  // cout << " paddle num = " << padnum << " ADC Pos hits = " << rawPosAdcHit.GetNPulses() << endl;
864  for (UInt_t thit=0; thit<rawPosAdcHit.GetNPulses(); ++thit) {
865  ((THcSignalHit*) frPosAdcPedRaw->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPedRaw());
866  ((THcSignalHit*) frPosAdcPed->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPed());
867 
868  ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseIntRaw(thit));
869  ((THcSignalHit*) frPosAdcPulseInt->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseInt(thit));
870 
871  ((THcSignalHit*) frPosAdcPulseAmpRaw->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseAmpRaw(thit));
872 
873  ((THcSignalHit*) frPosAdcPulseAmp->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseAmp(thit));
874 
875  ((THcSignalHit*) frPosAdcPulseTimeRaw->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseTimeRaw(thit));
876  ((THcSignalHit*) frPosAdcPulseTime->ConstructedAt(nrPosAdcHits))->Set(padnum, rawPosAdcHit.GetPulseTime(thit)+fAdcTdcOffset);
877 
878  if (rawPosAdcHit.GetPulseAmpRaw(thit) > 0) ((THcSignalHit*) frPosAdcErrorFlag->ConstructedAt(nrPosAdcHits))->Set(padnum, 0);
879  if (rawPosAdcHit.GetPulseAmpRaw(thit) <= 0) ((THcSignalHit*) frPosAdcErrorFlag->ConstructedAt(nrPosAdcHits))->Set(padnum, 1);
880 
881  ++nrPosAdcHits;
882  fTotNumAdcHits++;
883  fTotNumPosAdcHits++;
884  }
885  THcRawAdcHit& rawNegAdcHit = hit->GetRawAdcHitNeg();
886  if (rawNegAdcHit.GetNPulses()>0 && rawNegAdcHit.HasRefTime()) {
887  if (fNegAdcRefTime == kBig) {
888  fNegAdcRefTime=rawNegAdcHit.GetRefTime() ;
889  fNegAdcRefDiffTime=rawNegAdcHit.GetRefDiffTime() ;
890  }
891  if (fNegAdcRefTime != rawNegAdcHit.GetRefTime()) {
892  cout << "THcScintillatorPlane: " << GetName()<< " Neg ADC reftime problem at paddle num = " << padnum << " TDC pos hits = " << rawNegAdcHit.GetNPulses() << endl;
893  problem_flag=kTRUE;
894  }
895  }
896  // cout << " paddle num = " << padnum << " ADC Neg hits = " << rawNegAdcHit.GetNPulses() << endl;
897  for (UInt_t thit=0; thit<rawNegAdcHit.GetNPulses(); ++thit) {
898  ((THcSignalHit*) frNegAdcPedRaw->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPedRaw());
899  ((THcSignalHit*) frNegAdcPed->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPed());
900 
901  ((THcSignalHit*) frNegAdcPulseIntRaw->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseIntRaw(thit));
902  ((THcSignalHit*) frNegAdcPulseInt->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseInt(thit));
903 
904  ((THcSignalHit*) frNegAdcPulseAmpRaw->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseAmpRaw(thit));
905  ((THcSignalHit*) frNegAdcPulseAmp->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseAmp(thit));
906  ((THcSignalHit*) frNegAdcPulseTimeRaw->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseTimeRaw(thit));
907  ((THcSignalHit*) frNegAdcPulseTime->ConstructedAt(nrNegAdcHits))->Set(padnum, rawNegAdcHit.GetPulseTime(thit)+fAdcTdcOffset);
908 
909  if (rawNegAdcHit.GetPulseAmpRaw(thit) > 0) ((THcSignalHit*) frNegAdcErrorFlag->ConstructedAt(nrNegAdcHits))->Set(padnum, 0);
910  if (rawNegAdcHit.GetPulseAmpRaw(thit) <= 0) ((THcSignalHit*) frNegAdcErrorFlag->ConstructedAt(nrNegAdcHits))->Set(padnum, 1);
911 
912  ++nrNegAdcHits;
913  fTotNumAdcHits++;
914  fTotNumNegAdcHits++;
915  }
916 
917  // Need to be finding first hit in TDC range, not the first hit overall
918  if (hit->GetRawTdcHitPos().GetNHits() > 0)
919  ((THcSignalHit*) frPosTDCHits->ConstructedAt(nrPosTDCHits++))->Set(padnum, hit->GetRawTdcHitPos().GetTime()+fTdcOffset);
920  if (hit->GetRawTdcHitNeg().GetNHits() > 0)
921  ((THcSignalHit*) frNegTDCHits->ConstructedAt(nrNegTDCHits++))->Set(padnum, hit->GetRawTdcHitNeg().GetTime()+fTdcOffset);
922  // Should we make lists of offset corrected ADC Pulse times here too? For now
923  // the frNegAdcPulseTime frPosAdcPulseTime have that offset correction.
924  //
925  Bool_t badcraw_pos=kFALSE;
926  Bool_t badcraw_neg=kFALSE;
927  Double_t adcped_pos=-999;
928  Double_t adcped_neg=-999;
929  Int_t adcmult_pos=0;
930  Int_t adcmult_neg=0;
931  Int_t adchitused_pos=0;
932  Int_t adchitused_neg=0;
933  Double_t adcint_pos=-999;
934  Double_t adcint_neg=-999;
935  Double_t adcamp_pos=-kBig;
936  Double_t adcamp_neg=-kBig;
937  Double_t adctime_pos=kBig;
938  Double_t adctime_neg=kBig;
939  Double_t adctdcdifftime_pos=kBig;
940  Double_t adctdcdifftime_neg=kBig;
941  Double_t good_ielem_posadc = -1;
942  Double_t good_ielem_negadc = -1;
943  Bool_t btdcraw_pos=kFALSE;
944  Bool_t btdcraw_neg=kFALSE;
945  // Determine good tdc pos and neg times
946  Int_t tdc_pos=-999;
947  Int_t tdc_neg=-999;
948  Double_t good_ielem_postdc = -1;
949  Double_t good_ielem_negtdc = -1;
950  for(UInt_t thit=0; thit<hit->GetRawTdcHitPos().GetNHits(); thit++) {
951  tdc_pos = hit->GetRawTdcHitPos().GetTime(thit)+fTdcOffset;
952  if(tdc_pos >= fScinTdcMin && tdc_pos <= fScinTdcMax) {
953  btdcraw_pos = kTRUE;
954  good_ielem_postdc = thit;
955  break;
956  }
957  }
958  for(UInt_t thit=0; thit<hit->GetRawTdcHitNeg().GetNHits(); thit++) {
959  tdc_neg = hit->GetRawTdcHitNeg().GetTime(thit)+fTdcOffset;
960  if(tdc_neg >= fScinTdcMin && tdc_neg <= fScinTdcMax) {
961  btdcraw_neg = kTRUE;
962  good_ielem_negtdc = thit;
963  break;
964  }
965  }
966  //
967  if(fADCMode == kADCDynamicPedestal) {
968  Int_t good_ielem_negadc_test2=-1;
969  Int_t good_ielem_posadc_test2=-1;
970  //Loop Here over all hits per event for neg side of plane
971  if (good_ielem_negtdc != -1) {
972  Double_t max_adcamp_test=-1000.;
973  Double_t max_adctdcdiff_test=1000.;
974  for (UInt_t ielem=0;ielem<rawNegAdcHit.GetNPulses();ielem++) {
975  Double_t pulseAmp = rawNegAdcHit.GetPulseAmp(ielem);
976  Double_t pulseTime = rawNegAdcHit.GetPulseTime(ielem)+fAdcTdcOffset;
977  Double_t TdcAdcTimeDiff = tdc_neg*fScinTdcToTime-pulseTime;
978  if (rawNegAdcHit.GetPulseAmpRaw(ielem) <= 0)pulseAmp= 200.;
979  Bool_t pulseTimeCut =( TdcAdcTimeDiff > fHodoNegAdcTimeWindowMin[index]) && (TdcAdcTimeDiff < fHodoNegAdcTimeWindowMax[index]);
980  if (pulseTimeCut && pulseAmp>max_adcamp_test) {
981  good_ielem_negadc = ielem;
982  max_adcamp_test=pulseAmp;
983  }
984  if (abs(TdcAdcTimeDiff) < max_adctdcdiff_test) {
985  good_ielem_negadc_test2 = ielem;
986  max_adctdcdiff_test=abs(TdcAdcTimeDiff);
987  }
988  }
989  }
990 
991  //
992  if ( good_ielem_negadc == -1 && good_ielem_negadc_test2 != -1) good_ielem_negadc=good_ielem_negadc_test2;
993  if ( good_ielem_negadc == -1 && good_ielem_negadc_test2 == -1 && rawNegAdcHit.GetNPulses()>0) {
994  good_ielem_negadc=0;
995  }
996  //
997  if (good_ielem_negadc != -1 && good_ielem_negadc<rawNegAdcHit.GetNPulses()) {
998  adcped_neg = rawNegAdcHit.GetPed();
999  adcmult_neg = rawNegAdcHit.GetNPulses();
1000  adchitused_neg = good_ielem_negadc+1;
1001  adcint_neg = rawNegAdcHit.GetPulseInt(good_ielem_negadc);
1002  adcamp_neg = rawNegAdcHit.GetPulseAmp(good_ielem_negadc);
1003  if (rawNegAdcHit.GetPulseAmpRaw(good_ielem_negadc) <= 0) adcamp_neg= 200.;
1004  adctime_neg = rawNegAdcHit.GetPulseTime(good_ielem_negadc)+fAdcTdcOffset;
1005  badcraw_neg = kTRUE;
1006  adctdcdifftime_neg=tdc_neg*fScinTdcToTime-adctime_neg;
1007  }
1008  //
1009  //Loop Here over all hits per event for pos side of plane
1010  if (good_ielem_postdc != -1) {
1011  Double_t max_adcamp_test=-1000.;
1012  Double_t max_adctdcdiff_test=1000.;
1013  //
1014  for (UInt_t ielem=0;ielem<rawPosAdcHit.GetNPulses();ielem++) {
1015  Double_t pulseAmp = rawPosAdcHit.GetPulseAmp(ielem);
1016  Double_t pulseTime = rawPosAdcHit.GetPulseTime(ielem)+fAdcTdcOffset;
1017  Double_t TdcAdcTimeDiff = tdc_pos*fScinTdcToTime-pulseTime;
1018  Bool_t pulseTimeCut =( TdcAdcTimeDiff > fHodoPosAdcTimeWindowMin[index]) && (TdcAdcTimeDiff < fHodoPosAdcTimeWindowMax[index]);
1019  if (rawPosAdcHit.GetPulseAmpRaw(ielem) <= 0)pulseAmp= 200.;
1020  if (pulseTimeCut && pulseAmp>max_adcamp_test) {
1021  good_ielem_posadc = ielem;
1022  max_adcamp_test=pulseAmp;
1023  }
1024  if (abs(TdcAdcTimeDiff) < max_adctdcdiff_test) {
1025  good_ielem_posadc_test2 = ielem;
1026  max_adctdcdiff_test=abs(TdcAdcTimeDiff);
1027  }
1028  }
1029  } //
1030  if ( good_ielem_posadc == -1 && good_ielem_posadc_test2 != -1) good_ielem_posadc=good_ielem_posadc_test2;
1031  if ( good_ielem_posadc == -1 && good_ielem_posadc_test2 == -1 && rawPosAdcHit.GetNPulses()>0) good_ielem_posadc=0;
1032  if (good_ielem_posadc != -1 && good_ielem_posadc<rawPosAdcHit.GetNPulses()) {
1033  adcped_pos = rawPosAdcHit.GetPed();
1034  adcmult_pos = rawPosAdcHit.GetNPulses();
1035  adchitused_pos = good_ielem_posadc+1;
1036  adcint_pos = rawPosAdcHit.GetPulseInt(good_ielem_posadc);
1037  adcamp_pos = rawPosAdcHit.GetPulseAmp(good_ielem_posadc);
1038  if (rawPosAdcHit.GetPulseAmpRaw(good_ielem_posadc) <= 0) adcamp_pos= 200.;
1039  adctime_pos = rawPosAdcHit.GetPulseTime(good_ielem_posadc)+fAdcTdcOffset;
1040  badcraw_pos = kTRUE;
1041  adctdcdifftime_pos=tdc_pos*fScinTdcToTime-adctime_pos;
1042  //
1043  }
1044  } else if (fADCMode == kADCSampleIntegral) {
1045  adcint_pos = hit->GetRawAdcHitPos().GetSampleIntRaw() - fPosPed[index];
1046  adcint_neg = hit->GetRawAdcHitNeg().GetSampleIntRaw() - fNegPed[index];
1047  badcraw_pos = badcraw_neg = kTRUE;
1048  } else if (fADCMode == kADCSampIntDynPed) {
1049  adcint_pos = hit->GetRawAdcHitPos().GetSampleInt();
1050  adcint_neg = hit->GetRawAdcHitNeg().GetSampleInt();
1051  badcraw_pos = badcraw_neg = kTRUE;
1052  } else {
1053  adcint_pos = hit->GetRawAdcHitPos().GetPulseIntRaw()-fPosPed[index];
1054  adcint_neg = hit->GetRawAdcHitNeg().GetPulseIntRaw()-fNegPed[index];
1055  badcraw_pos = badcraw_neg = kTRUE;
1056  }
1057  if (adcint_pos >= fADCDiagCut) {
1058  ((THcSignalHit*) frPosADCHits->ConstructedAt(nrPosADCHits))->Set(padnum, adcint_pos);
1059  Double_t samplesum=hit->GetRawAdcHitPos().GetSampleIntRaw();
1060  Double_t pedestal=hit->GetRawAdcHitPos().GetPedRaw();
1061  ((THcSignalHit*) frPosADCSums->ConstructedAt(nrPosADCHits))->Set(padnum, samplesum);
1062  ((THcSignalHit*) frPosADCPeds->ConstructedAt(nrPosADCHits++))->Set(padnum, pedestal);
1063  }
1064  if (adcint_neg >= fADCDiagCut) {
1065  ((THcSignalHit*) frNegADCHits->ConstructedAt(nrNegADCHits))->Set(padnum, adcint_neg);
1066  Double_t samplesum=hit->GetRawAdcHitNeg().GetSampleIntRaw();
1067  Double_t pedestal=hit->GetRawAdcHitNeg().GetPedRaw();
1068  ((THcSignalHit*) frNegADCSums->ConstructedAt(nrNegADCHits))->Set(padnum, samplesum);
1069  ((THcSignalHit*) frNegADCPeds->ConstructedAt(nrNegADCHits++))->Set(padnum, pedestal);
1070  }
1071  //
1072  if((btdcraw_pos && badcraw_pos) || (btdcraw_neg && badcraw_neg )) {
1073  if (good_ielem_posadc != -1) {
1074  //good adc multiplicity
1075  fTotNumGoodPosAdcHits++;
1076  fTotNumGoodAdcHits++;
1077  //good adc occupancy
1078  fNumGoodPosAdcHits.at(padnum-1) = padnum;
1079  fGoodPosAdcPed.at(padnum-1) = adcped_pos;
1080  fGoodPosAdcMult.at(padnum-1) = adcmult_pos;
1081  fGoodPosAdcHitUsed.at(padnum-1) = adchitused_pos;
1082  fGoodPosAdcPulseInt.at(padnum-1) = adcint_pos;
1083  fGoodPosAdcPulseAmp.at(padnum-1) = adcamp_pos;
1084  fGoodPosAdcPulseTime.at(padnum-1) = adctime_pos;
1085  fGoodPosAdcTdcDiffTime.at(padnum-1) = adctdcdifftime_pos;
1086  }
1087  if (good_ielem_negadc != -1) {
1088  //good adc multiplicity
1089  fTotNumGoodNegAdcHits++;
1090  fTotNumGoodAdcHits++;
1091  //good adc occupancy
1092  fNumGoodNegAdcHits.at(padnum-1) = padnum;
1093  fGoodNegAdcPed.at(padnum-1) = adcped_neg;
1094  fGoodNegAdcMult.at(padnum-1) = adcmult_neg;
1095  fGoodNegAdcHitUsed.at(padnum-1) = adchitused_neg;
1096  fGoodNegAdcPulseInt.at(padnum-1) = adcint_neg;
1097  fGoodNegAdcPulseAmp.at(padnum-1) = adcamp_neg;
1098  fGoodNegAdcPulseTime.at(padnum-1) = adctime_neg;
1099  fGoodNegAdcTdcDiffTime.at(padnum-1) = adctdcdifftime_neg;
1100  }
1101 
1102  //DEFINE THE "GOOD +TDC Multiplicities and Occupancies"
1103  if (good_ielem_postdc != -1) {
1104  fTotNumGoodPosTdcHits++;
1105  fTotNumGoodTdcHits++;
1106  //good tdc occupancy
1107  fNumGoodPosTdcHits.at(padnum-1) = padnum;
1108  }
1109 
1110  //DEFINE THE "GOOD -TDC Multiplicities and Occupancies"
1111  if (good_ielem_negtdc != -1) {
1112  fTotNumGoodNegTdcHits++;
1113  fTotNumGoodTdcHits++;
1114  //good tdc occupancy
1115  fNumGoodNegTdcHits.at(padnum-1) = padnum;
1116  }
1117  new( (*fHodoHits)[fNScinHits]) THcHodoHit(tdc_pos, tdc_neg,
1118  adcint_pos, adcint_neg,
1119  hit->fCounter, this);
1120  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCpeak(adcamp_pos);
1121  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCpeak(adcamp_neg);
1122  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCtime(adctime_pos);
1123  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCtime(adctime_neg);
1124 
1125  //CALCULATE TIME-WALK CORRECTIONS HERE!!!!
1126  //Define GoodTdcUnCorrTime
1127  if(btdcraw_pos&&badcraw_pos) {
1128  fGoodPosTdcTimeUnCorr.at(padnum-1) = tdc_pos*fScinTdcToTime;
1129  //tw_corr_pos = fHodoPos_c1[padnum-1]/pow(adcamp_pos/fTdc_Thrs,fHodoPos_c2[padnum-1]) - fHodoPos_c1[padnum-1]/pow(200./fTdc_Thrs, fHodoPos_c2[padnum-1]);
1130 
1131  tw_corr_pos = 1./pow(adcamp_pos/fTdc_Thrs,fHodoPos_c2[padnum-1]) - 1./pow(200./fTdc_Thrs, fHodoPos_c2[padnum-1]);
1132 
1133  fGoodPosTdcTimeWalkCorr.at(padnum-1) = tdc_pos*fScinTdcToTime -tw_corr_pos;
1134  }
1135  if(btdcraw_neg&&badcraw_neg) {
1136  fGoodNegTdcTimeUnCorr.at(padnum-1) = tdc_neg*fScinTdcToTime;
1137 
1138  //tw_corr_neg = fHodoNeg_c1[padnum-1]/pow(adcamp_neg/fTdc_Thrs,fHodoNeg_c2[padnum-1]) - fHodoNeg_c1[padnum-1]/pow(200./fTdc_Thrs, fHodoNeg_c2[padnum-1]);
1139 
1140  tw_corr_neg = 1./pow(adcamp_neg/fTdc_Thrs,fHodoNeg_c2[padnum-1]) - 1./pow(200./fTdc_Thrs, fHodoNeg_c2[padnum-1]);
1141 
1142  fGoodNegTdcTimeWalkCorr.at(padnum-1) = tdc_neg*fScinTdcToTime -tw_corr_neg;
1143 
1144  }
1145 
1146 
1147  // Do corrections if valid TDC on both ends of bar
1148  if( (btdcraw_pos && btdcraw_neg) && (badcraw_pos && badcraw_neg) ) {
1149  // Do the pulse height correction to the time. (Position dependent corrections later)
1150  Double_t adc_timec_pos= adctime_pos;
1151  Double_t adc_timec_neg= adctime_neg;
1152  Double_t timec_pos, timec_neg;
1153  if(fTofUsingInvAdc) {
1154  timec_pos = tdc_pos*fScinTdcToTime
1155  - fHodoPosInvAdcOffset[index]
1156  - fHodoPosInvAdcAdc[index]/TMath::Sqrt(TMath::Max(20.0*.020,adcint_pos));
1157  timec_neg = tdc_neg*fScinTdcToTime
1158  - fHodoNegInvAdcOffset[index]
1159  - fHodoNegInvAdcAdc[index]/TMath::Sqrt(TMath::Max(20.0*.020,adcint_neg));
1160  } else { // FADC style
1161  timec_pos = tdc_pos*fScinTdcToTime -tw_corr_pos + fHodo_LCoeff[index];
1162  timec_neg = tdc_neg*fScinTdcToTime -tw_corr_neg- 2*fHodoCableFit[index] + fHodo_LCoeff[index];
1163  adc_timec_pos = adc_timec_pos -tw_corr_pos + fHodo_LCoeff[index];
1164  adc_timec_neg = adc_timec_neg -tw_corr_neg- 2*fHodoCableFit[index] + fHodo_LCoeff[index];
1165  }
1166 
1167  Double_t TWCorrDiff = fGoodNegTdcTimeWalkCorr.at(padnum-1) - 2*fHodoCableFit[index] - fGoodPosTdcTimeWalkCorr.at(padnum-1);
1168 
1169  Double_t fHitDistCorr = 0.5*TWCorrDiff*fHodoVelFit[index];
1170 
1171 
1172 
1173  fGoodDiffDistTrack.at(index) = fHitDistCorr;
1174 
1175  Double_t vellight=fHodoVelLight[index]; //read from hodo_cuts.param, where it is set fixed to 15.0
1176 
1177  Double_t dist_from_center=0.5*(timec_neg-timec_pos)*vellight;
1178  Double_t scint_center=0.5*(fPosLeft+fPosRight);
1179  Double_t hit_position=scint_center+dist_from_center;
1180  hit_position=TMath::Min(hit_position,fPosLeft);
1181  hit_position=TMath::Max(hit_position,fPosRight);
1182  Double_t scin_corrected_time, postime, negtime;
1183  Double_t adc_postime=adc_timec_pos;
1184  Double_t adc_negtime=adc_timec_neg;
1185  if(fTofUsingInvAdc) {
1186  timec_pos -= (fPosLeft-hit_position)/
1187  fHodoPosInvAdcLinear[index];
1188  timec_neg -= (hit_position-fPosRight)/
1189  fHodoNegInvAdcLinear[index];
1190  scin_corrected_time = 0.5*(timec_pos+timec_neg);
1191  if (fCosmicFlag) {
1192  postime = timec_pos + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1193  negtime = timec_neg + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1194  } else {
1195  postime = timec_pos - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1196  negtime = timec_neg - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1197  }
1198  } else {
1199  scin_corrected_time = 0.5*(timec_neg+timec_pos);
1200  timec_pos= scin_corrected_time;
1201  timec_neg= scin_corrected_time;
1202  Double_t adc_time_corrected = 0.5*(adc_timec_pos+adc_timec_neg);
1203  if (fCosmicFlag) {
1204  postime = timec_pos + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1205  negtime = timec_neg + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1206  adc_postime = adc_time_corrected + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1207  adc_negtime = adc_time_corrected + (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1208  } else {
1209  postime = timec_pos - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1210  negtime = timec_neg - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1211  adc_postime = adc_time_corrected - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1212  adc_negtime = adc_time_corrected - (fZpos+(index%2)*fDzpos)/(29.979*fBetaNominal);
1213  }
1214  }
1215  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPaddleCenter(fPosCenter[index]);
1216  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetCorrectedTimes(timec_pos,timec_neg,
1217  postime, negtime,
1218  scin_corrected_time);
1219  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCpeak(adcamp_pos);
1220  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCpeak(adcamp_neg);
1221  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCCorrtime(adc_postime);
1222  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCCorrtime(adc_negtime);
1223  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetCalcPosition(fHitDistCorr); //
1224 
1225  fGoodPosTdcTimeCorr.at(padnum-1) = timec_pos;
1226  fGoodNegTdcTimeCorr.at(padnum-1) = timec_neg;
1227  fGoodPosTdcTimeTOFCorr.at(padnum-1) = postime;
1228  fGoodNegTdcTimeTOFCorr.at(padnum-1) = negtime;
1229  } else {
1230  Double_t timec_pos,timec_neg;
1231  timec_pos=tdc_pos;
1232  timec_neg=tdc_neg;
1233  if(btdcraw_pos&& badcraw_pos) {
1234  if(fTofUsingInvAdc) {
1235  timec_pos = tdc_pos*fScinTdcToTime
1236  - fHodoPosInvAdcOffset[index]
1237  - fHodoPosInvAdcAdc[index]/TMath::Sqrt(TMath::Max(20.0*.020,adcint_pos));
1238  } else { // FADC style
1239  timec_pos = tdc_pos*fScinTdcToTime -tw_corr_pos + fHodo_LCoeff[index];
1240  }
1241  }
1242  if(btdcraw_neg && badcraw_neg) {
1243  if(fTofUsingInvAdc) {
1244  timec_neg = tdc_neg*fScinTdcToTime
1245  - fHodoNegInvAdcOffset[index]
1246  - fHodoNegInvAdcAdc[index]/TMath::Sqrt(TMath::Max(20.0*.020,adcint_neg));
1247  } else { // FADC style
1248  timec_neg = tdc_neg*fScinTdcToTime -tw_corr_neg- 2*fHodoCableFit[index] + fHodo_LCoeff[index];
1249  }
1250  }
1251  Double_t adc_neg=0.,adc_pos=0.;
1252  if (badcraw_neg) adc_neg=adcamp_neg;
1253  if (badcraw_pos) adc_pos=adcamp_pos;
1254  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPaddleCenter(fPosCenter[index]);
1255  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetCorrectedTimes(timec_pos,timec_neg);
1256  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCpeak(adc_neg); // needed for new TWCOrr
1257  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCpeak(adc_pos); // needed for new TWCOrr
1258  if (badcraw_neg) {
1259  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCtime(adctime_neg);
1260  } else {
1261  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetNegADCtime(-999.);
1262  }
1263  if (badcraw_pos) {
1264  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCtime(adctime_pos);
1265  } else {
1266  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetPosADCtime(-999.);
1267  }
1268  ((THcHodoHit*) fHodoHits->At(fNScinHits))->SetCalcPosition(kBig); //
1269  fGoodPosTdcTimeCorr.at(padnum-1) = timec_pos;
1270  fGoodNegTdcTimeCorr.at(padnum-1) = timec_neg;
1271  fGoodPosTdcTimeTOFCorr.at(padnum-1) = kBig;
1272  fGoodNegTdcTimeTOFCorr.at(padnum-1) = kBig;
1273  }
1274  // if ( ((THcHodoHit*) fHodoHits->At(fNScinHits))->GetPosTOFCorrectedTime() != ((THcHodoHit*) fHodoHits->At(fNScinHits))->GetPosTOFCorrectedTime()) cout << " ihit = " << ihit << " scinhit = " << fNScinHits << " plane = " << fPlaneNum << " padnum = " << padnum << " " << tdc_pos<< " "<< tdc_neg<< " " << ((THcHodoHit*) fHodoHits->At(fNScinHits))->GetPosTOFCorrectedTime() << endl;
1275  fNScinHits++; // One or more good time counter
1276  //
1277  }
1278  ihit++; // Raw hit counter
1279  }
1280  // cout << "THcScintillatorPlane: ihit = " << ihit << endl;
1281  if (problem_flag) {
1282  cout << "THcScintillatorPlane::ProcessHits " << fPlaneNum << " " << nexthit << "/" << nrawhits << endl;
1283 cout << " Ref problem end *******" << endl;
1284  }
1285  return(ihit);
1286 }
1287 
1288 //_____________________________________________________________________________
1290 {
1295  Int_t nrawhits = rawhits->GetLast()+1;
1296  // cout << "THcScintillatorPlane::AcculatePedestals " << fPlaneNum << " " << nexthit << "/" << nrawhits << endl;
1297 
1298  Int_t ihit = nexthit;
1299  while(ihit < nrawhits) {
1300  THcRawHodoHit* hit = (THcRawHodoHit *) rawhits->At(ihit);
1301  if(hit->fPlane > fPlaneNum) {
1302  break;
1303  }
1304  Int_t element = hit->fCounter - 1; // Should check if in range
1305  Int_t adcpos = hit->GetRawAdcHitPos().GetPulseIntRaw();
1306  Int_t adcneg = hit->GetRawAdcHitNeg().GetPulseIntRaw();
1307 
1308  if(adcpos <= fPosPedLimit[element]) {
1309  fPosPedSum[element] += adcpos;
1310  fPosPedSum2[element] += adcpos*adcpos;
1311  fPosPedCount[element]++;
1312  if(fPosPedCount[element] == fMinPeds/5) {
1313  fPosPedLimit[element] = 100 + fPosPedSum[element]/fPosPedCount[element];
1314  }
1315  }
1316  if(adcneg <= fNegPedLimit[element]) {
1317  fNegPedSum[element] += adcneg;
1318  fNegPedSum2[element] += adcneg*adcneg;
1319  fNegPedCount[element]++;
1320  if(fNegPedCount[element] == fMinPeds/5) {
1321  fNegPedLimit[element] = 100 + fNegPedSum[element]/fNegPedCount[element];
1322  }
1323  }
1324  ihit++;
1325  }
1326 
1327  fNPedestalEvents++;
1328 
1329  return(ihit);
1330 }
1331 
1332 //_____________________________________________________________________________
1334 {
1340  for(UInt_t i=0; i<fNelem;i++) {
1341 
1342  // Positive tubes
1343  fPosPed[i] = ((Double_t) fPosPedSum[i]) / TMath::Max(1, fPosPedCount[i]);
1344  fPosThresh[i] = fPosPed[i] + 15;
1345 
1346  // Negative tubes
1347  fNegPed[i] = ((Double_t) fNegPedSum[i]) / TMath::Max(1, fNegPedCount[i]);
1348  fNegThresh[i] = fNegPed[i] + 15;
1349 
1350  // cout <<"Pedestals "<< i+1 << " " << fPosPed[i] << " " << fNegPed[i] << endl;
1351  }
1352  // cout << " " << endl;
1353 
1354 }
1355 
1356 //_____________________________________________________________________________
1358 {
1364  fNPedestalEvents = 0;
1365  fMinPeds = 500; // In engine, this is set in parameter file
1366  fPosPedSum = new Int_t [fNelem];
1367  fPosPedSum2 = new Int_t [fNelem];
1368  fPosPedLimit = new Int_t [fNelem];
1369  fPosPedCount = new Int_t [fNelem];
1370  fNegPedSum = new Int_t [fNelem];
1371  fNegPedSum2 = new Int_t [fNelem];
1372  fNegPedLimit = new Int_t [fNelem];
1373  fNegPedCount = new Int_t [fNelem];
1374 
1375  fPosPed = new Double_t [fNelem];
1376  fNegPed = new Double_t [fNelem];
1377  fPosThresh = new Double_t [fNelem];
1378  fNegThresh = new Double_t [fNelem];
1379  for(UInt_t i=0;i<fNelem;i++) {
1380  fPosPedSum[i] = 0;
1381  fPosPedSum2[i] = 0;
1382  fPosPedLimit[i] = 1000; // In engine, this are set in parameter file
1383  fPosPedCount[i] = 0;
1384  fNegPedSum[i] = 0;
1385  fNegPedSum2[i] = 0;
1386  fNegPedLimit[i] = 1000; // In engine, this are set in parameter file
1387  fNegPedCount[i] = 0;
1388  }
1389 }
1390 //____________________________________________________________________________
Int_t GetTimeRaw(UInt_t iHit=0) const
Gets raw TDC time. In channels.
std::string GetName(const std::string &scope_name)
A single plane of scintillators.
Double_t GetSampleInt() const
Gets pedestal subtracted integral of samples. In channels.
Double_t GetTofTolerance() const
Definition: THcHodoscope.h:68
Double_t GetHodoNegInvAdcAdc(Int_t iii) const
Definition: THcHodoscope.h:83
Double_t GetHodoNegAdcTimeWindowMax(Int_t iii) const
Definition: THcHodoscope.h:86
Int_t GetScinIndex(Int_t nPlane, Int_t nPaddle)
Int_t GetLast() const
const char Option_t
Class representing a single hit for the Hodoscopes.
Definition: THcHodoHit.h:16
Double_t GetHodoPosTimeOffset(Int_t iii) const
Definition: THcHodoscope.h:75
virtual Int_t ReadDatabase(const TDatime &date)
Double_t GetHodoPosAdcTimeWindowMax(Int_t iii) const
Definition: THcHodoscope.h:84
Double_t GetHodoNegTimeOffset(Int_t iii) const
Definition: THcHodoscope.h:76
Int_t GetPulseTimeRaw(UInt_t iPulse=0) const
Gets raw pulse time. In subsamples.
Short_t Min(Short_t a, Short_t b)
int Int_t
bool Bool_t
Int_t GetTdcOffset(Int_t ip) const
Definition: THcHodoscope.h:116
Double_t GetHodoCableFit(Int_t iii) const
Definition: THcHodoscope.h:91
Double_t GetHodoNeg_c1(Int_t iii) const
Definition: THcHodoscope.h:96
STL namespace.
UInt_t GetNHits() const
Gets the number of set hits.
Class representing a single raw TDC hit.
Definition: THcRawTdcHit.h:7
Double_t GetAdcTdcOffset(Int_t ip) const
Definition: THcHodoscope.h:117
Int_t fCounter
Definition: THcRawHit.h:55
Double_t GetHodoVelFit(Int_t iii) const
Definition: THcHodoscope.h:90
Double_t GetStartTimeSlop() const
Definition: THcHodoscope.h:102
Double_t GetHodoPos_c2(Int_t iii) const
Definition: THcHodoscope.h:97
Double_t GetTdcMin() const
Definition: THcHodoscope.h:66
double pow(double, double)
Double_t GetHodoLCoeff(Int_t iii) const
Definition: THcHodoscope.h:92
Double_t GetHodoNegAdcTimeWindowMin(Int_t iii) const
Definition: THcHodoscope.h:87
Double_t GetPed() const
Gets sample pedestal. In channels.
virtual Int_t AccumulatePedestals(TClonesArray *rawhits, Int_t nexthit)
THcRawAdcHit & GetRawAdcHitNeg()
Int_t GetRefDiffTime() const
THcRawTdcHit & GetRawTdcHitNeg()
THcRawAdcHit & GetRawAdcHitPos()
Double_t GetHodoNeg_c2(Int_t iii) const
Definition: THcHodoscope.h:98
Int_t GetPulseAmpRaw(UInt_t iPulse=0) const
Gets raw pulse amplitude. In channels.
Double_t GetHodoNegSigma(Int_t iii) const
Definition: THcHodoscope.h:124
Double_t GetPulseAmp(UInt_t iPulse=0) const
Gets pedestal subtracted pulse amplitude. In channels.
Double_t GetPulseInt(UInt_t iPulse=0) const
Gets pedestal subtracted pulse integral. In channels.
Double_t GetHodoNegPhcCoeff(Int_t iii) const
Definition: THcHodoscope.h:72
virtual Int_t CoarseProcess(TClonesArray &tracks)
Bool_t HasRefTime() const
Queries whether reference time has been set.
Double_t GetTdcMax() const
Definition: THcHodoscope.h:67
Double_t GetTdcToTime() const
Definition: THcHodoscope.h:69
unsigned int UInt_t
Int_t GetSampleIntRaw() const
Gets raw integral of sTimeFacamples. In channels.
char * Form(const char *fmt,...)
Double_t GetHodoPos_c1(Int_t iii) const
Definition: THcHodoscope.h:95
virtual void Clear(Option_t *opt="")
tuple list
Definition: SConscript.py:9
const Bool_t kFALSE
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
virtual Int_t DefineVariables(EMode mode=kDefine)
Double_t GetHodoPosInvAdcLinear(Int_t iii) const
Definition: THcHodoscope.h:80
Int_t GetRefTime() const
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
double Double_t
Bool_t HasRefTime() const
Double_t GetHodoNegInvAdcOffset(Int_t iii) const
Definition: THcHodoscope.h:79
Int_t GetTime(UInt_t iHit=0) const
Gets TDC time. In channels.
Double_t GetHodoSlop(Int_t ip)
Definition: THcHodoscope.h:113
virtual Int_t ProcessHits(TClonesArray *rawhits, Int_t nexthit)
UInt_t GetNPulses() const
Gets number of set pulses.
virtual EStatus Init(const TDatime &run_time)
Double_t GetHodoVelLight(Int_t iii) const
Definition: THcHodoscope.h:77
Int_t fPlane
Definition: THcRawHit.h:54
Double_t GetHodoPosInvAdcAdc(Int_t iii) const
Definition: THcHodoscope.h:82
Double_t GetBetaNominal() const
Definition: THcHodoscope.h:70
THcRawTdcHit & GetRawTdcHitPos()
Double_t GetStartTimeCenter() const
Definition: THcHodoscope.h:101
ClassImp(THcScintillatorPlane) THcScintillatorPlane
Int_t GetRefTime() const
Gets reference time. In channels.
Int_t GetPedRaw() const
Gets raw signal pedestal. In channels.
Short_t Max(Short_t a, Short_t b)
Double_t GetHodoNegMinPh(Int_t iii) const
Definition: THcHodoscope.h:74
Double_t GetHodoPosSigma(Int_t iii) const
Definition: THcHodoscope.h:123
Double_t GetPulseTime(UInt_t iPulse=0) const
R__EXTERN class THcParmList * gHcParms
Definition: THcGlobals.h:11
Double_t GetHodoPosAdcTimeWindowMin(Int_t iii) const
Definition: THcHodoscope.h:85
Double_t GetHodoPosPhcCoeff(Int_t iii) const
Definition: THcHodoscope.h:71
Class representing a single raw hit for a hodoscope paddle.
Definition: THcRawHodoHit.h:9
Double_t Sqrt(Double_t x)
TObject * At(Int_t idx) const
static const Double_t kBig
Definition: THcFormula.cxx:31
Int_t GetRefDiffTime() const
Double_t GetHodoNegInvAdcLinear(Int_t iii) const
Definition: THcHodoscope.h:81
const Bool_t kTRUE
Generic hodoscope consisting of multiple planes with multiple paddles with phototubes on both ends...
Definition: THcHodoscope.h:37
Double_t GetHodoPosInvAdcOffset(Int_t iii) const
Definition: THcHodoscope.h:78
Double_t GetHodoPosMinPh(Int_t iii) const
Definition: THcHodoscope.h:73
Class representing a single raw ADC hit.
Definition: THcRawAdcHit.h:7
virtual Int_t FineProcess(TClonesArray &tracks)
virtual Int_t Decode(const THaEvData &)
Double_t GetTDCThrs() const
Definition: THcHodoscope.h:99