Hall C ROOT/C++ Analyzer (hcana)
THcHodoscope.cxx
Go to the documentation of this file.
1 
12 #include "THcSignalHit.h"
13 #include "THcScintPlaneCluster.h"
14 #include "THcShower.h"
15 #include "THcCherenkov.h"
16 #include "THcHallCSpectrometer.h"
17 #include "THcHitList.h"
18 #include "THcRawShowerHit.h"
19 #include "TClass.h"
20 #include "math.h"
21 #include "THaSubDetector.h"
22 
23 #include "THcHodoscope.h"
24 #include "THaEvData.h"
25 #include "THaDetMap.h"
26 #include "THcDetectorMap.h"
27 #include "THaGlobals.h"
28 #include "THaCutList.h"
29 #include "THcGlobals.h"
30 #include "THcParmList.h"
31 #include "VarDef.h"
32 #include "VarType.h"
33 #include "THaTrack.h"
34 #include "TClonesArray.h"
35 #include "TMath.h"
36 
37 #include "THaTrackProj.h"
38 #include <vector>
39 
40 #include <cstring>
41 #include <cstdio>
42 #include <cstdlib>
43 #include <iostream>
44 #include <iomanip>
45 #include <fstream>
46 #include <cassert>
47 
48 using namespace std;
49 using std::vector;
50 
51 //_____________________________________________________________________________
52 THcHodoscope::THcHodoscope( const char* name, const char* description,
53  THaApparatus* apparatus ) :
54  THaNonTrackingDetector(name,description,apparatus)
55 {
56  // Constructor
57 
58  //fTrackProj = new TClonesArray( "THaTrackProj", 5 );
59  // Construct the planes
60  fNPlanes = 0; // No planes until we make them
61  fStartTime=-1e5;
63 }
64 
65 //_____________________________________________________________________________
67  THaNonTrackingDetector()
68 {
69  // Constructor
70 }
71 
72 //_____________________________________________________________________________
73 void THcHodoscope::Setup(const char* name, const char* description)
74 {
84  if( IsZombie()) return;
85 
86  // fDebug = 1; // Keep this at one while we're working on the code
87 
88  char prefix[2];
89 
90  prefix[0]=tolower(GetApparatus()->GetName()[0]);
91  prefix[1]='\0';
92 
93  TString temp(prefix[0]);
94  fSHMS=kFALSE;
95  if (temp == "p" ) fSHMS=kTRUE;
96  TString histname=temp+"_timehist";
97  hTime = new TH1F(histname,"",400,0,200);
98  // cout << " fSHMS = " << fSHMS << endl;
99  string planenamelist;
100  DBRequest listextra[]={
101  {"hodo_num_planes", &fNPlanes, kInt},
102  {"hodo_plane_names",&planenamelist, kString},
103  {"hodo_tdcrefcut", &fTDC_RefTimeCut, kInt, 0, 1},
104  {"hodo_adcrefcut", &fADC_RefTimeCut, kInt, 0, 1},
105  {0}
106  };
107  //fNPlanes = 4; // Default if not defined
108  fTDC_RefTimeCut = 0; // Minimum allowed reference times
109  fADC_RefTimeCut = 0;
110  gHcParms->LoadParmValues((DBRequest*)&listextra,prefix);
111 
112  cout << "Plane Name List : " << planenamelist << endl;
113 
114  vector<string> plane_names = Podd::vsplit(planenamelist);
115  // Plane names
116  if(plane_names.size() != (UInt_t) fNPlanes) {
117  cout << "ERROR: Number of planes " << fNPlanes << " doesn't agree with number of plane names " << plane_names.size() << endl;
118  // Should quit. Is there an official way to quit?
119  }
120  fPlaneNames = new char* [fNPlanes];
121  for(Int_t i=0;i<fNPlanes;i++) {
122  fPlaneNames[i] = new char[plane_names[i].length()+1];
123  strcpy(fPlaneNames[i], plane_names[i].c_str());
124  }
125 
126  // Probably shouldn't assume that description is defined
127  char* desc = new char[strlen(description)+100];
129  for(Int_t i=0;i < fNPlanes;i++) {
130  strcpy(desc, description);
131  strcat(desc, " Plane ");
132  strcat(desc, fPlaneNames[i]);
133  fPlanes[i] = new THcScintillatorPlane(fPlaneNames[i], desc, i+1, this); // Number planes starting from zero!!
134  cout << "Created Scintillator Plane " << fPlaneNames[i] << ", " << desc << endl;
135  }
136 
137  // Save the nominal particle mass
138  THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
139  fPartMass = app->GetParticleMass();
141 
142 
143 
144  if (fSHMS) {
145  fCherenkov = dynamic_cast<THcCherenkov*>(app->GetDetector("hgcer"));
146  } else {
147  fCherenkov = dynamic_cast<THcCherenkov*>(app->GetDetector("cer"));
148  }
149 
150  delete [] desc;
151 }
152 
153 //_____________________________________________________________________________
154 THaAnalysisObject::EStatus THcHodoscope::Init( const TDatime& date )
155 {
156  // cout << "In THcHodoscope::Init()" << endl;
157  Setup(GetName(), GetTitle());
158 
159  char EngineDID[] = "xSCIN";
160  EngineDID[0] = toupper(GetApparatus()->GetName()[0]);
161  if( gHcDetectorMap->FillMap(fDetMap, EngineDID) < 0 ) {
162  static const char* const here = "Init()";
163  Error( Here(here), "Error filling detectormap for %s.", EngineDID );
164  return kInitError;
165  }
166 
167  // Should probably put this in ReadDatabase as we will know the
168  // maximum number of hits after setting up the detector map
169  // But it needs to happen before the sub detectors are initialized
170  // so that they can get the pointer to the hitlist.
171  cout << " Hodo tdc ref time cut = " << fTDC_RefTimeCut << " " << fADC_RefTimeCut << endl;
172 
173  InitHitList(fDetMap, "THcRawHodoHit", fDetMap->GetTotNumChan()+1,
175 
176  EStatus status;
177  // This triggers call of ReadDatabase and DefineVariables
178  if( (status = THaNonTrackingDetector::Init( date )) )
179  return fStatus=status;
180 
181  for(Int_t ip=0;ip<fNPlanes;ip++) {
182  if((status = fPlanes[ip]->Init( date ))) {
183  return fStatus=status;
184  }
185  }
186 
187  fNScinHits = new Int_t [fNPlanes];
189  fNPlaneTime = new Int_t [fNPlanes];
191 
192  // Double_t fHitCnt4 = 0., fHitCnt3 = 0.;
193 
194  // Int_t m = 0;
195  // fScinHit = new Double_t*[fNPlanes];
196  // for ( m = 0; m < fNPlanes; m++ ){
197  // fScinHit[m] = new Double_t[fNPaddle[0]];
198  // }
199 
200  for (int ip=0; ip<fNPlanes; ++ip) {
201  fScinHitPaddle.push_back(std::vector<Int_t>(fNPaddle[ip], 0));
202  }
203 
204  fPresentP = 0;
205  THaVar* vpresent = gHaVars->Find(Form("%s.present",GetApparatus()->GetName()));
206  if(vpresent) {
207  fPresentP = (Bool_t *) vpresent->GetValuePointer();
208  }
209 
210  return fStatus = kOK;
211 }
212 //_____________________________________________________________________________
214 {
221  // static const char* const here = "ReadDatabase()";
222  char prefix[2];
223  char parname[100];
224 
225  // Determine which spectrometer in order to construct
226  // the parameter names (e.g. hscin_1x_nr vs. sscin_1x_nr)
227 
228  prefix[0]=tolower(GetApparatus()->GetName()[0]);
229  //
230  prefix[1]='\0';
231  strcpy(parname,prefix);
232  strcat(parname,"scin_");
233  //
234  //
235  // Int_t plen=strlen(parname);
236  // cout << " readdatabse hodo fnplanes = " << fNPlanes << endl;
237 
238  CreateMissReportParms(Form("%sscin",prefix));
239 
240  fBetaNoTrk = 0.;
241  fBetaNoTrkChiSq = 0.;
242  fNPaddle = new UInt_t [fNPlanes];
243  fFPTime = new Double_t [fNPlanes];
246 
247  prefix[0]=tolower(GetApparatus()->GetName()[0]);
248  //
249  prefix[1]='\0';
250 
251  for(Int_t i=0;i<fNPlanes;i++) {
252 
253  DBRequest list[]={
254  {Form("scin_%s_nr",fPlaneNames[i]), &fNPaddle[i], kInt},
255  {0}
256  };
257 
258 
259  gHcParms->LoadParmValues((DBRequest*)&list,prefix);
260  }
261 
262  // GN added
263  // reading variables from *hodo.param
266  for (Int_t i=1;i<fNPlanes;i++) {
268  fTotHodScin+=(fNPaddle[i]);
269  }
270  // need this for "padded arrays" i.e. 4x16 lists of parameters (GN)
272  if (fDebug>=1) cout <<"fMaxScinPerPlane = "<<fMaxScinPerPlane<<" fMaxHodoScin = "<<fMaxHodoScin<<endl;
273 
291 
292  //New Time-Walk Calibration Parameters
302 
303  fNHodoscopes = 2;
304  fxLoScin = new Int_t [fNHodoscopes];
305  fxHiScin = new Int_t [fNHodoscopes];
306  fyLoScin = new Int_t [fNHodoscopes];
307  fyHiScin = new Int_t [fNHodoscopes];
308  fHodoSlop = new Double_t [fNPlanes];
309  fTdcOffset = new Int_t [fNPlanes];
315 
316 
317  for(Int_t ip=0;ip<fNPlanes;ip++) { // Set a large default window
318  fTdcOffset[ip] = 0 ;
319  fAdcTdcOffset[ip] = 0.0 ;
320  }
321 
322  DBRequest list[]={
323  {"cosmicflag", &fCosmicFlag, kInt, 0, 1},
324  {"NumPlanesBetaCalc", &fNumPlanesBetaCalc, kInt, 0, 1},
325  {"start_time_center", &fStartTimeCenter, kDouble},
326  {"start_time_slop", &fStartTimeSlop, kDouble},
327  {"scin_tdc_to_time", &fScinTdcToTime, kDouble},
328  {"scin_tdc_min", &fScinTdcMin, kDouble},
329  {"scin_tdc_max", &fScinTdcMax, kDouble},
330  {"tof_tolerance", &fTofTolerance, kDouble, 0, 1},
331  {"pathlength_central", &fPathLengthCentral, kDouble},
332  {"hodo_pos_sigma", &fHodoPosSigma[0], kDouble, fMaxHodoScin, 1},
333  {"hodo_neg_sigma", &fHodoNegSigma[0], kDouble, fMaxHodoScin, 1},
334  {"hodo_pos_ped_limit", &fHodoPosPedLimit[0], kInt, fMaxHodoScin, 1},
335  {"hodo_neg_ped_limit", &fHodoNegPedLimit[0], kInt, fMaxHodoScin, 1},
336  {"tofusinginvadc", &fTofUsingInvAdc, kInt, 0, 1},
337  {"xloscin", &fxLoScin[0], kInt, (UInt_t) fNHodoscopes},
338  {"xhiscin", &fxHiScin[0], kInt, (UInt_t) fNHodoscopes},
339  {"yloscin", &fyLoScin[0], kInt, (UInt_t) fNHodoscopes},
340  {"yhiscin", &fyHiScin[0], kInt, (UInt_t) fNHodoscopes},
341  {"track_eff_test_num_scin_planes", &fTrackEffTestNScinPlanes, kInt},
342  {"trackeff_scint_ydiff_max", &trackeff_scint_ydiff_max, kDouble, 0, 1},
343  {"trackeff_scint_xdiff_max", &trackeff_scint_xdiff_max, kDouble, 0, 1},
344  {"cer_npe", &fNCerNPE, kDouble, 0, 1},
345  {"normalized_energy_tot", &fNormETot, kDouble, 0, 1},
346  {"hodo_slop", fHodoSlop, kDouble, (UInt_t) fNPlanes},
347  {"debugprintscinraw", &fdebugprintscinraw, kInt, 0,1},
348  {"hodo_tdc_offset", fTdcOffset, kInt, (UInt_t) fNPlanes, 1},
349  {"hodo_adc_tdc_offset", fAdcTdcOffset, kDouble, (UInt_t) fNPlanes, 1},
350  {"hodo_PosAdcTimeWindowMin", fHodoPosAdcTimeWindowMin, kDouble, (UInt_t) fMaxHodoScin, 1},
351  {"hodo_PosAdcTimeWindowMax", fHodoPosAdcTimeWindowMax, kDouble, (UInt_t) fMaxHodoScin, 1},
352  {"hodo_NegAdcTimeWindowMin", fHodoNegAdcTimeWindowMin, kDouble, (UInt_t) fMaxHodoScin, 1},
353  {"hodo_NegAdcTimeWindowMax", fHodoNegAdcTimeWindowMax, kDouble, (UInt_t) fMaxHodoScin, 1},
354  {"dumptof", &fDumpTOF, kInt, 0, 1},
355  {"TOFCalib_shtrk_lo", &fTOFCalib_shtrk_lo, kDouble, 0, 1},
356  {"TOFCalib_shtrk_hi", &fTOFCalib_shtrk_hi, kDouble, 0, 1},
357  {"TOFCalib_cer_lo", &fTOFCalib_cer_lo, kDouble, 0, 1},
358  {"TOFCalib_beta_lo", &fTOFCalib_beta_lo, kDouble, 0, 1},
359  {"TOFCalib_beta_hi", &fTOFCalib_beta_hi, kDouble, 0, 1},
360  {"dumptof_filename", &fTOFDumpFile, kString, 0, 1},
361  {"TrackBetaIncludeSinglePmtHits", &fTrackBetaIncludeSinglePmtHits, kInt, 0, 1},
362  {0}
363  };
364 
365  // Defaults if not defined in parameter file
366  fTrackBetaIncludeSinglePmtHits=0; // do not use paddles with only one hit in the TRack Beta calculation set ==1 to include
369  for(UInt_t ip=0;ip<fMaxHodoScin;ip++) {
370  fHodoPosAdcTimeWindowMin[ip] = -1000.;
371  fHodoPosAdcTimeWindowMax[ip] = 1000.;
372  fHodoNegAdcTimeWindowMin[ip] = -1000.;
373  fHodoNegAdcTimeWindowMax[ip] = 1000.;
374 
375  fHodoPosPedLimit[ip] = 0.0;
376  fHodoNegPedLimit[ip] = 0.0;
377  fHodoPosSigma[ip] = 0.2;
378  fHodoNegSigma[ip] = 0.2;
379  }
386  fDumpTOF = 0;
387  fTOFDumpFile="";
388  fTofUsingInvAdc = 1;
389  fTofTolerance = 3.0;
390  fNCerNPE = 2.0;
391  fNormETot = 0.7;
392  fCosmicFlag=0;
394  // Gets added to each reference time corrected raw TDC value
395  // to make sure valid range is all positive.
396 
397  gHcParms->LoadParmValues((DBRequest*)&list,prefix);
398  if (fCosmicFlag==1) cout << "Setup for cosmics in TOF"<< endl;
399  // cout << " cosmic flag = " << fCosmicFlag << endl;
400  if(fDumpTOF) {
401  fDumpOut.open(fTOFDumpFile.c_str());
402  if(fDumpOut.is_open()) {
403  //fDumpOut << "Hodoscope Time of Flight calibration data" << endl;
404  } else {
405  fDumpTOF = 0;
406  cout << "WARNING: Unable to open TOF Dump file " << fTOFDumpFile << endl;
407  cout << "Data for TOF calibration not being written." << endl;
408  }
409  }
410 
411  // cout << " x1 lo = " << fxLoScin[0]
412  // << " x2 lo = " << fxLoScin[1]
413  // << " x1 hi = " << fxHiScin[0]
414  // << " x2 hi = " << fxHiScin[1]
415  // << endl;
416 
417  // cout << " y1 lo = " << fyLoScin[0]
418  // << " y2 lo = " << fyLoScin[1]
419  // << " y1 hi = " << fyHiScin[0]
420  // << " y2 hi = " << fyHiScin[1]
421  // << endl;
422 
423  // cout << "Hdososcope planes hits for trigger = " << fTrackEffTestNScinPlanes
424  // << " normalized energy min = " << fNormETot
425  // << " number of photo electrons = " << fNCerNPE
426  // << endl;
427 
428  //Set scin Velocity/Cable to default
429  for (UInt_t i=0; i<fMaxHodoScin; i++) {
430  fHodoVelLight[i] = 15.0;
431  }
432 
433  if (fTofUsingInvAdc) {
434  DBRequest list2[]={
435  {"hodo_vel_light", &fHodoVelLight[0], kDouble, fMaxHodoScin, 1},
436  {"hodo_pos_invadc_offset",&fHodoPosInvAdcOffset[0],kDouble,fMaxHodoScin},
437  {"hodo_neg_invadc_offset",&fHodoNegInvAdcOffset[0],kDouble,fMaxHodoScin},
438  {"hodo_pos_invadc_linear",&fHodoPosInvAdcLinear[0],kDouble,fMaxHodoScin},
439  {"hodo_neg_invadc_linear",&fHodoNegInvAdcLinear[0],kDouble,fMaxHodoScin},
440  {"hodo_pos_invadc_adc",&fHodoPosInvAdcAdc[0],kDouble,fMaxHodoScin},
441  {"hodo_neg_invadc_adc",&fHodoNegInvAdcAdc[0],kDouble,fMaxHodoScin},
442  {0}
443  };
444 
445  gHcParms->LoadParmValues((DBRequest*)&list2,prefix);
446  };
447  /* if (!fTofUsingInvAdc) {
448  DBRequest list3[]={
449  {"hodo_vel_light", &fHodoVelLight[0], kDouble, fMaxHodoScin},
450  {"hodo_pos_minph", &fHodoPosMinPh[0], kDouble, fMaxHodoScin},
451  {"hodo_neg_minph", &fHodoNegMinPh[0], kDouble, fMaxHodoScin},
452  {"hodo_pos_phc_coeff", &fHodoPosPhcCoeff[0], kDouble, fMaxHodoScin},
453  {"hodo_neg_phc_coeff", &fHodoNegPhcCoeff[0], kDouble, fMaxHodoScin},
454  {"hodo_pos_time_offset", &fHodoPosTimeOffset[0], kDouble, fMaxHodoScin},
455  {"hodo_neg_time_offset", &fHodoNegTimeOffset[0], kDouble, fMaxHodoScin},
456  {0}
457  };
458 
459 
460  gHcParms->LoadParmValues((DBRequest*)&list3,prefix);
461  */
462  DBRequest list4[]={
463  {"hodo_velFit", &fHodoVelFit[0], kDouble, fMaxHodoScin, 1},
464  {"hodo_cableFit", &fHodoCableFit[0], kDouble, fMaxHodoScin, 1},
465  {"hodo_LCoeff", &fHodo_LCoeff[0], kDouble, fMaxHodoScin, 1},
466  {"c1_Pos", &fHodoPos_c1[0], kDouble, fMaxHodoScin, 1},
467  {"c1_Neg", &fHodoNeg_c1[0], kDouble, fMaxHodoScin, 1},
468  {"c2_Pos", &fHodoPos_c2[0], kDouble, fMaxHodoScin, 1},
469  {"c2_Neg", &fHodoNeg_c2[0], kDouble, fMaxHodoScin, 1},
470  {"TDC_threshold", &fTdc_Thrs, kDouble, 0, 1},
471  {"hodo_PosSigma", &fHodoSigmaPos[0], kDouble, fMaxHodoScin, 1},
472  {"hodo_NegSigma", &fHodoSigmaNeg[0], kDouble, fMaxHodoScin, 1},
473  {0}
474  };
475 
476  fTdc_Thrs = 1.0;
477  //Set Default Values if NOT defined in param file
478  for (UInt_t i=0; i<fMaxHodoScin; i++)
479  {
480 
481  //Turn OFF Time-Walk Correction if param file NOT found
482  fHodoPos_c1[i] = 0.0;
483  fHodoPos_c2[i] = 0.0;
484  fHodoNeg_c1[i] = 0.0;
485  fHodoNeg_c2[i] = 0.0;
486  }
487  for (UInt_t i=0; i<fMaxHodoScin; i++)
488  {
489  //Set scin Velocity/Cable to default
490  fHodoCableFit[i] = 0.0;
491  fHodoVelFit[i] = 15.0;
492  //set time coeff between paddles to default
493  fHodo_LCoeff[i] = 0.0;
494 
495  }
496 
497  gHcParms->LoadParmValues((DBRequest*)&list4,prefix);
498 
499  if (fDebug >=1) {
500  cout <<"******* Testing Hodoscope Parameter Reading ***\n";
501  cout<<"StarTimeCenter = "<<fStartTimeCenter<<endl;
502  cout<<"StartTimeSlop = "<<fStartTimeSlop<<endl;
503  cout <<"ScintTdcToTime = "<<fScinTdcToTime<<endl;
504  cout <<"TdcMin = "<<fScinTdcMin<<" TdcMax = "<<fScinTdcMax<<endl;
505  cout <<"TofTolerance = "<<fTofTolerance<<endl;
506  cout <<"*** VelLight ***\n";
507  for (Int_t i1=0;i1<fNPlanes;i1++) {
508  cout<<"Plane "<<i1<<endl;
509  for (UInt_t i2=0;i2<fMaxScinPerPlane;i2++) {
510  cout<<fHodoVelLight[GetScinIndex(i1,i2)]<<" ";
511  }
512  cout <<endl;
513  }
514  cout <<endl<<endl;
515  // check fHodoPosPhcCoeff
516  /*
517  cout <<"fHodoPosPhcCoeff = ";
518  for (int i1=0;i1<fMaxHodoScin;i1++) {
519  cout<<this->GetHodoPosPhcCoeff(i1)<<" ";
520  }
521  cout<<endl;
522  */
523  }
524  //
525  if ((fTofTolerance > 0.5) && (fTofTolerance < 10000.)) {
526  cout << "USING "<<fTofTolerance<<" NSEC WINDOW FOR FP NO_TRACK CALCULATIONS.\n";
527  }
528  else {
529  fTofTolerance= 3.0;
530  cout << "*** USING DEFAULT 3 NSEC WINDOW FOR FP NO_TRACK CALCULATIONS!! ***\n";
531  }
532  //
533  fRatio_xpfp_to_xfp=0.00;
534  TString SHMS="p";
535  TString HMS="h";
536  TString test=prefix[0];
537  if (test==SHMS ) fRatio_xpfp_to_xfp=0.0018; // SHMS
538  if (test == HMS ) fRatio_xpfp_to_xfp=0.0011; // HMS
539  cout << " fRatio_xpfp_to_xfp= " << fRatio_xpfp_to_xfp << endl;
540  //
541  fIsInit = true;
542  return kOK;
543 }
544 
545 //_____________________________________________________________________________
547 {
551  // cout << "THcHodoscope::DefineVariables called " << GetName() << endl;
552  if( mode == kDefine && fIsSetup ) return kOK;
553  fIsSetup = ( mode == kDefine );
554 
555  // Register variables in global list
556 
557  RVarDef vars[] = {
558  // Move these into THcHallCSpectrometer using track fTracks
559  {"beta", "Beta including track info", "fBeta"},
560  {"betanotrack", "Beta from scintillator hits", "fBetaNoTrk"},
561  {"betachisqnotrack", "Chi square of beta from scintillator hits", "fBetaNoTrkChiSq"},
562  {"fpHitsTime", "Time at focal plane from all hits", "fFPTimeAll"},
563  {"starttime", "Hodoscope Start Time", "fStartTime"},
564  {"goodstarttime", "Hodoscope Good Start Time (logical flag)", "fGoodStartTime"},
565  {"goodscinhit", "Hit in fid area", "fGoodScinHits"},
566  {"adctdc_offset"," ","fOffsetTime"},
567  {"TimeHist_StartTime_Sigma", "", "fTimeHist_StartTime_Sigma"},
568  {"TimeHist_StartTime_Peak", "", "fTimeHist_StartTime_Peak"},
569  {"TimeHist_StartTime_NumPeaks", "", "fTimeHist_StartTime_NumPeaks"},
570  {"TimeHist_StartTime_Hits", "", "fTimeHist_StartTime_Hits"},
571  {"TimeHist_FpTime_Sigma", "", "fTimeHist_FpTime_Sigma"},
572  {"TimeHist_FpTime_Peak", "", "fTimeHist_FpTime_Peak"},
573  {"TimeHist_FpTime_NumPeaks", "", "fTimeHist_FpTime_NumPeaks"},
574  {"TimeHist_FpTime_Hits", "", "fTimeHist_FpTime_Hits"},
575  { 0 }
576  };
577  return DefineVarsFromList( vars, mode );
578  // return kOK;
579 }
580 
581 //_____________________________________________________________________________
583 {
584  // Destructor. Remove variables from global list.
585 
586  delete [] fFPTime;
587  delete [] fPlaneCenter;
588  delete [] fPlaneSpacing;
589 
590  if( fIsSetup )
591  RemoveVariables();
592  if( fIsInit )
593  DeleteArrays();
594 
595  for( int i = 0; i < fNPlanes; ++i ) {
596  delete fPlanes[i];
597  delete [] fPlaneNames[i];
598  }
599  delete [] fPlanes;
600  delete [] fPlaneNames;
601 }
602 
603 //_____________________________________________________________________________
605 {
606  // Delete member arrays. Used by destructor.
607  // Int_t k;
608  // for( k = 0; k < fNPlanes; k++){
609  // delete [] fScinHit[k];
610  // }
611  // delete [] fScinHit;
612 
613  delete [] fxLoScin; fxLoScin = NULL;
614  delete [] fxHiScin; fxHiScin = NULL;
615  delete [] fyLoScin; fyLoScin = NULL;
616  delete [] fyHiScin; fyHiScin = NULL;
617  delete [] fHodoSlop; fHodoSlop = NULL;
618 
619  delete [] fNPaddle; fNPaddle = NULL;
620  delete [] fHodoVelLight; fHodoVelLight = NULL;
621  delete [] fHodoPosSigma; fHodoPosSigma = NULL;
622  delete [] fHodoNegSigma; fHodoNegSigma = NULL;
623  delete [] fHodoPosMinPh; fHodoPosMinPh = NULL;
624  delete [] fHodoNegMinPh; fHodoNegMinPh = NULL;
625  delete [] fHodoPosPhcCoeff; fHodoPosPhcCoeff = NULL;
626  delete [] fHodoNegPhcCoeff; fHodoNegPhcCoeff = NULL;
627  delete [] fHodoPosTimeOffset; fHodoPosTimeOffset = NULL;
628  delete [] fHodoNegTimeOffset; fHodoNegTimeOffset = NULL;
629  delete [] fHodoPosPedLimit; fHodoPosPedLimit = NULL;
630  delete [] fHodoNegPedLimit; fHodoNegPedLimit = NULL;
631  delete [] fHodoPosInvAdcOffset; fHodoPosInvAdcOffset = NULL;
632  delete [] fHodoNegInvAdcOffset; fHodoNegInvAdcOffset = NULL;
633  delete [] fHodoPosInvAdcLinear; fHodoPosInvAdcLinear = NULL;
634  delete [] fHodoNegInvAdcLinear; fHodoNegInvAdcLinear = NULL;
635  delete [] fHodoPosInvAdcAdc; fHodoPosInvAdcAdc = NULL;
636  delete [] fHodoNegInvAdcAdc; fHodoNegInvAdcAdc = NULL;
637  delete [] fGoodPlaneTime; fGoodPlaneTime = NULL;
638  delete [] fNPlaneTime; fNPlaneTime = NULL;
639  delete [] fSumPlaneTime; fSumPlaneTime = NULL;
640  delete [] fNScinHits; fNScinHits = NULL;
641  delete [] fTdcOffset; fTdcOffset = NULL;
642  delete [] fAdcTdcOffset; fAdcTdcOffset = NULL;
647 
648  delete [] fHodoVelFit; fHodoVelFit = NULL;
649  delete [] fHodoCableFit; fHodoCableFit = NULL;
650  delete [] fHodo_LCoeff; fHodo_LCoeff = NULL;
651  delete [] fHodoPos_c1; fHodoPos_c1 = NULL;
652  delete [] fHodoNeg_c1; fHodoNeg_c1 = NULL;
653  delete [] fHodoPos_c2; fHodoPos_c2 = NULL;
654  delete [] fHodoNeg_c2; fHodoNeg_c2 = NULL;
655  delete [] fHodoSigmaPos; fHodoSigmaPos = NULL;
656  delete [] fHodoSigmaNeg; fHodoSigmaNeg = NULL;
657 }
658 
659 //_____________________________________________________________________________
661 {
675 
676  fBeta = 0.0;
677  fBetaNoTrk = 0.0;
678  fBetaNoTrkChiSq = 0.0;
679  fStartTime = -1000.;
680  fADCStartTime = -1000.;
681  fOffsetTime = kBig;
682  fFPTimeAll= -1000.;
684  fGoodScinHits = 0;
685 
686  if( *opt != 'I' ) {
687  for(Int_t ip=0;ip<fNPlanes;ip++) {
688  fPlanes[ip]->Clear();
689  fFPTime[ip]=0.;
690  fPlaneCenter[ip]=0.;
691  fPlaneSpacing[ip]=0.;
692  for(UInt_t iPaddle=0;iPaddle<fNPaddle[ip]; ++iPaddle) {
693  fScinHitPaddle[ip][iPaddle]=0;
694  }
695  }
696  }
697  fdEdX.clear();
698  fNScinHit.clear();
699  fNClust.clear();
700  fClustSize.clear();
701  fClustPos.clear();
702  fNCluster.clear();
703  fClusterSize.clear();
704  fClusterXPos.clear();
705  fClusterYPos.clear();
706  fThreeScin.clear();
707  fGoodScinHitsX.clear();
708  fGoodFlags.clear();
709 }
710 
711 //_____________________________________________________________________________
713 {
726  // Get the Hall C style hitlist (fRawHitList) for this event
727  Bool_t present = kTRUE; // Suppress reference time warnings
728  if(fPresentP) { // if this spectrometer not part of trigger
729  present = *fPresentP;
730  }
731  fNHits = DecodeToHitList(evdata, !present);
732  fEventNum = evdata.GetEvNum();
733  //
734  // GN: print event number so we can cross-check with engine
735  // if (evdata.GetEvNum()>1000)
736  // cout <<"\nhcana_event " << evdata.GetEvNum()<<endl;
737 
738  fCheckEvent = evdata.GetEvNum();
739  fEventType = evdata.GetEvType();
740 
741  if(gHaCuts->Result("Pedestal_event")) {
742  Int_t nexthit = 0;
743  for(Int_t ip=0;ip<fNPlanes;ip++) {
744 
745  nexthit = fPlanes[ip]->AccumulatePedestals(fRawHitList, nexthit);
746  }
747  fAnalyzePedestals = 1; // Analyze pedestals first normal events
748  return(0);
749  }
750  if(fAnalyzePedestals) {
751  for(Int_t ip=0;ip<fNPlanes;ip++) {
752 
754  }
755  fAnalyzePedestals = 0; // Don't analyze pedestals next event
756  }
757 
758  // Let each plane get its hits
759  Int_t nexthit = 0;
760  //THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
761  // cout << " event number = " << fEventNum << " Evtyp = " << fEventType<< " spec = " << app->GetName() << endl;
762  fNfptimes=0;
763  Int_t thits = 0;
764  for(Int_t ip=0;ip<fNPlanes;ip++) {
765 
766  fPlaneCenter[ip] = fPlanes[ip]->GetPosCenter(0) + fPlanes[ip]->GetPosOffset();
767  fPlaneSpacing[ip] = fPlanes[ip]->GetSpacing();
768 
769  // nexthit = fPlanes[ip]->ProcessHits(fRawHitList, nexthit);
770  // GN: select only events that have reasonable TDC values to start with
771  // as per the Engine h_strip_scin.f
772  nexthit = fPlanes[ip]->ProcessHits(fRawHitList,nexthit);
773  thits+=fPlanes[ip]->GetNScinHits();
774  }
775  //
776  //
777  fStartTime=-1000;
778  if (thits>0 ) EstimateFocalPlaneTime();
779 
780  if (fdebugprintscinraw == 1) {
781  for(UInt_t ihit = 0; ihit < fNRawHits ; ihit++) {
782 // THcRawHodoHit* hit = (THcRawHodoHit *) fRawHitList->At(ihit);
783 // cout << ihit << " : " << hit->fPlane << ":" << hit->fCounter << " : "
784 // << hit->fADC_pos << " " << hit->fADC_neg << " " << hit->fTDC_pos
785 // << " " << hit->fTDC_neg << endl;
786  }
787  cout << endl;
788  }
789 
790 
791  return fNHits;
792 }
793 //_____________________________________________________________________________
795 {
796  Double_t time_peak=-1000;
797  Int_t NBinsX=hTime->GetNbinsX();
798  Int_t hTimeScanRange = 10.; // Integrate over HtimeScanRange
799  vector<Double_t> hpeakCent;
800  vector<Double_t> hpeakNum;
801  vector<Double_t> hpeakRMS;
802  vector<Double_t> hpeakFlag;
803  vector<Int_t> hpeakBin;
804  Double_t MinimumNum=2.;
805  Double_t test_peakmax=0.;
806  Bool_t scanning_for_local_peak=kFALSE;
807  Double_t save_mean=0,save_rms=0,save_num=0;
808  Int_t save_bin;
809  Bool_t new_peak=kFALSE;
810  Bool_t replace_peak=kFALSE;
811  UInt_t best_peak_index=0;
812  Int_t best_peak_num=-1;
813  Double_t best_peak_diff=1000;
814  Int_t nfound=0;
815  for (Int_t nb=1;nb<NBinsX-hTimeScanRange;nb++) {
816  hTime->GetXaxis()->SetRange(nb,nb+hTimeScanRange);
817  Double_t test_int = hTime->Integral();
818  if (scanning_for_local_peak) {
819  if ( test_int <= test_peakmax) {
820  Int_t ps=hpeakCent.size();
821  replace_peak=kFALSE;
822  new_peak=kFALSE;
823  if (ps==0) new_peak=kTRUE;
824  if (ps!=0 && nb==hpeakBin[ps]+1 && save_mean!=hpeakCent[ps-1]) new_peak=kFALSE;
825  if (ps!=0 && save_num > hpeakNum[ps-1] && abs(save_mean-hpeakCent[ps-1])<5) replace_peak=kTRUE;
826  if (ps!=0 && nb!=hpeakBin[ps]+1 && save_num > MinimumNum && abs(save_mean-hpeakCent[ps-1])>=5) new_peak=kTRUE;
827  if (new_peak) {
828  hpeakCent.push_back(save_mean);
829  hpeakRMS.push_back(save_rms);
830  hpeakNum.push_back(save_num);
831  hpeakBin.push_back(save_bin);
832  hpeakFlag.push_back(1);
833  }
834  if (replace_peak) {
835  hpeakCent[ps-1]=save_mean;
836  hpeakRMS[ps-1]=save_rms;
837  hpeakNum[ps-1]=save_num;
838  hpeakBin[ps-1]=save_bin;
839  hpeakFlag[ps-1]=1;
840  }
841  scanning_for_local_peak = kFALSE;
842  test_peakmax = 0;
843  best_peak_index=-1;
844  best_peak_num=5;
845  nfound=0;
846  for (UInt_t np=0;np<hpeakNum.size();np++) {
847  hpeakFlag[np]=-1;
848  if ( hpeakNum[np] > 5 && (hpeakNum[np]>= best_peak_num || abs(hpeakNum[np] - best_peak_num)<= 4) ) {
849  if (nfound==0 || (hpeakNum[np]== best_peak_num || abs(hpeakNum[np] - best_peak_num)<= 4) ) {
850  hpeakFlag[np]=1;
851  if (nfound==0 ) best_peak_num =hpeakNum[np] ;
852  if (nfound==0 ) best_peak_index= np;
853  nfound++;
854  } else {
855  for (UInt_t nt=0;nt<np;nt++) {hpeakFlag[nt]=-1;}
856  nfound=1;
857  best_peak_num =hpeakNum[np] ;
858  best_peak_index= np;
859  }
860  }
861  }
862  if (nfound>1) {
863  best_peak_diff=1000;
864  for (UInt_t np=0;np<hpeakNum.size();np++) {
865  if (hpeakFlag[np]==1 && abs(hpeakCent[np]-fStartTimeCenter)<best_peak_diff) {
866  best_peak_diff = abs(hpeakCent[np]-fStartTimeCenter);
867  best_peak_index= np;
868  }
869  }}
870  if (nfound==0) {
871  best_peak_diff=1000;
872  for (UInt_t np=0;np<hpeakNum.size();np++) {
873  if (abs(hpeakCent[np]-fStartTimeCenter)<best_peak_diff) {
874  best_peak_diff = abs(hpeakCent[np]-fStartTimeCenter);
875  best_peak_index= np;
876  }
877  }}
878 
879  } else {
880  test_peakmax = test_int;
881  save_mean=hTime->GetMean();
882  save_rms=hTime->GetRMS();
883  save_num=hTime->Integral();
884  save_bin=nb;
885  }
886  } else {
887  if ( test_int > MinimumNum) {
888  test_peakmax = test_int;
889  scanning_for_local_peak = kTRUE;
890  save_mean=hTime->GetMean();
891  save_rms=hTime->GetRMS();
892  save_num=hTime->Integral();
893  save_bin=nb;
894  }
895  }
896  }
897  //
898  if (hpeakNum.size() >0 && best_peak_index<hpeakNum.size() ) {
899  time_peak= hpeakCent[best_peak_index];
900  if (FillFlag==1) {
901  fTimeHist_StartTime_NumPeaks=hpeakNum.size() ;
902  fTimeHist_StartTime_Peak= time_peak;
903  fTimeHist_StartTime_Sigma= hpeakRMS[best_peak_index] ;
904  fTimeHist_StartTime_Hits= hpeakNum[best_peak_index];
905  }
906  if (FillFlag==2) {
907  fTimeHist_FpTime_NumPeaks=hpeakNum.size() ;
908  fTimeHist_FpTime_Peak= time_peak;
909  fTimeHist_FpTime_Sigma= hpeakRMS[best_peak_index] ;
910  fTimeHist_FpTime_Hits= hpeakNum[best_peak_index];
911  }
912  }
913  //
914  return time_peak;
915 }
916 
917 //_____________________________________________________________________________
919 {
929  Int_t ihit=0;
930  Int_t nscinhits=0; // Total # hits with at least one good tdc
931  hTime->Reset();
932  //
933  //
934  for(Int_t ip=0;ip<fNPlanes;ip++) {
935  Int_t nphits=fPlanes[ip]->GetNScinHits();
936  nscinhits += nphits;
937  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
938  for(Int_t i=0;i<nphits;i++) {
939  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(i);
940  if(hit->GetHasCorrectedTimes()) {
941  Double_t postime=hit->GetPosTOFCorrectedTime();
942  Double_t negtime=hit->GetNegTOFCorrectedTime();
943  hTime->Fill(postime);
944  hTime->Fill(negtime);
945  }
946  }
947  }
948  //
949  Double_t TimePeak=DetermineTimePeak(1);
950  hTime->Reset();
951  //
952  Double_t AdcTdcDiffTimeSum=0;
953  Double_t NAdcTdcDiffTimeSum=0;
954  //
955  for(Int_t ip=0;ip<fNPlanes;ip++) {
956  Int_t nphits=fPlanes[ip]->GetNScinHits();
957  nscinhits += nphits;
958  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
959  for(Int_t i=0;i<nphits;i++) {
960  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(i);
961  if(hit->GetHasCorrectedTimes()) {
962  NAdcTdcDiffTimeSum++;
963  AdcTdcDiffTimeSum+=(hit->GetPosADCtime()-hit->GetPosTDC()*fScinTdcToTime);
964  NAdcTdcDiffTimeSum++;
965  AdcTdcDiffTimeSum+=(hit->GetNegADCtime()-hit->GetNegTDC()*fScinTdcToTime);
966  Double_t postime=hit->GetPosADCCorrtime();
967  Double_t negtime=hit->GetNegADCCorrtime();
968  hTime->Fill(postime);
969  hTime->Fill(negtime);
970  }
971  }
972  }
973  if (NAdcTdcDiffTimeSum>0) AdcTdcDiffTimeSum=AdcTdcDiffTimeSum/NAdcTdcDiffTimeSum;
974  //
975  Double_t AdcTimePeak=DetermineTimePeak(3);
976  //
977  ihit = 0;
978  Double_t fpTimeSum = 0.0;
979  Double_t adcfpTimeSum = 0.0;
980  Double_t adcNfptimes=0;
981  fNfptimes=0;
982  Int_t Ngood_hits_plane=0;
983  Int_t Ngood_adchits_plane=0;
984  Double_t Plane_fptime_sum=0.0;
985  Double_t Plane_adcfptime_sum=0.0;
986  Bool_t goodplanetime[fNPlanes];
987  Bool_t twogoodtimes[nscinhits];
988  Int_t NumPlanesGoodHit=0;
989  Int_t NumPlanesGoodAdcHit=0;
990  if (TimePeak>0) {
991  for(Int_t ip=0;ip<fNumPlanesBetaCalc;ip++) {
992  goodplanetime[ip] = kFALSE;
993  Int_t nphits=fPlanes[ip]->GetNScinHits();
994  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
995  Ngood_hits_plane=0;
996  Plane_fptime_sum=0.0;
997  for(Int_t i=0;i<nphits;i++) {
998  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(i);
999  twogoodtimes[ihit] = kFALSE;
1000  if(hit->GetHasCorrectedTimes()) {
1001  Double_t postime=hit->GetPosTOFCorrectedTime();
1002  Double_t negtime=hit->GetNegTOFCorrectedTime();
1003  Double_t adcpostime=hit->GetPosADCCorrtime();
1004  Double_t adcnegtime=hit->GetNegADCCorrtime();
1005  if ((postime>(TimePeak-fTofTolerance)) && (postime<(TimePeak+fTofTolerance)) &&
1006  (negtime>(TimePeak-fTofTolerance)) && (negtime<(TimePeak+fTofTolerance)) ) {
1007  hit->SetTwoGoodTimes(kTRUE);
1008  twogoodtimes[ihit] = kTRUE; // Both tubes fired
1009  Int_t index=hit->GetPaddleNumber()-1; //
1010  Double_t fptime;
1011  if(fCosmicFlag==1) {
1012  fptime = hit->GetScinCorrectedTime()
1013  + (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
1014  / (29.979 * fBetaNominal);
1015  }else{
1016  fptime = hit->GetScinCorrectedTime()
1017  - (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
1018  / (29.979 * fBetaNominal);
1019  }
1020  Ngood_hits_plane++;
1021  Plane_fptime_sum+=fptime;
1022  fpTimeSum += fptime;
1023  fNfptimes++;
1024  goodplanetime[ip] = kTRUE;
1025  } else {
1026  hit->SetTwoGoodTimes(kFALSE);
1027  }
1028  //
1029  if ((adcpostime>(AdcTimePeak-fTofTolerance)) && (adcpostime<(AdcTimePeak+fTofTolerance)) &&
1030  (adcnegtime>(AdcTimePeak-fTofTolerance)) && (adcnegtime<(AdcTimePeak+fTofTolerance)) ) {
1031  Int_t index=hit->GetPaddleNumber()-1; //
1032  Double_t fptime;
1033  if(fCosmicFlag==1) {
1034  fptime = hit->GetScinCorrectedTime()
1035  + (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
1036  / (29.979 * fBetaNominal);
1037  }else{
1038  fptime = hit->GetScinCorrectedTime()
1039  - (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
1040  / (29.979 * fBetaNominal);
1041  }
1042  Ngood_adchits_plane++;
1043  Plane_adcfptime_sum+=fptime;
1044  adcfpTimeSum += fptime;
1045  adcNfptimes++;
1046  }
1047  //
1048  }
1049  ihit++;
1050  }
1051  if (Ngood_hits_plane>0) NumPlanesGoodHit++;
1052  if (Ngood_adchits_plane>0) NumPlanesGoodAdcHit++;
1053  if (Ngood_hits_plane>0) fPlanes[ip]->SetFpTime(Plane_fptime_sum/float(Ngood_hits_plane));
1054  fPlanes[ip]->SetNGoodHits(Ngood_hits_plane);
1055  }
1056  } // if TimePeak>0
1057  //
1058  if(NumPlanesGoodHit>=3) {
1059  fStartTime = fpTimeSum/fNfptimes;
1062  fOffsetTime=0;
1063  if(NumPlanesGoodAdcHit>=3) {
1064  fADCStartTime = adcfpTimeSum/adcNfptimes-fStartTime;
1065  }
1066  fOffsetTime =AdcTdcDiffTimeSum;
1067  } else {
1072  fOffsetTime=AdcTdcDiffTimeSum;
1073  }
1074  //
1075  //
1076  hTime->Reset();
1077  //
1078  if(fGoodStartTime && (goodplanetime[0]||goodplanetime[1]) &&(goodplanetime[2]||goodplanetime[3])) {
1079 
1080  Double_t sumW = 0.;
1081  Double_t sumT = 0.;
1082  Double_t sumZ = 0.;
1083  Double_t sumZZ = 0.;
1084  Double_t sumTZ = 0.;
1085  Int_t ihhit = 0;
1086 
1087  for(Int_t ip=0;ip<fNumPlanesBetaCalc;ip++) {
1088  Int_t nphits=fPlanes[ip]->GetNScinHits();
1089  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
1090 
1091  for(Int_t i=0;i<nphits;i++) {
1092  Int_t index=((THcHodoHit*)hodoHits->At(i))->GetPaddleNumber()-1;
1093 
1094  if(twogoodtimes[ihhit]){
1095 
1096  Double_t sigma = 0.0;
1097  if(fTofUsingInvAdc)
1098  {
1099  sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoPosSigma[GetScinIndex(ip,index)],2) +
1100  TMath::Power( fHodoNegSigma[GetScinIndex(ip,index)],2) ) );
1101  }
1102  else{
1103  sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoSigmaPos[GetScinIndex(ip,index)],2) +
1104  TMath::Power( fHodoSigmaNeg[GetScinIndex(ip,index)],2) ) );
1105  }
1106 
1107  Double_t scinWeight = 1 / TMath::Power(sigma,2);
1108  Double_t zPosition = fPlanes[ip]->GetZpos() + (index%2)*fPlanes[ip]->GetDzpos();
1109  // cout << "hit = " << ihhit + 1 << " zpos = " << zPosition << " sigma = " << sigma << endl;
1110  //cout << "fHodoSigma+ = " << fHodoSigmaPos[GetScinIndex(ip,index)] << endl;
1111  sumW += scinWeight;
1112  sumT += scinWeight * ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime();
1113  sumZ += scinWeight * zPosition;
1114  sumZZ += scinWeight * ( zPosition * zPosition );
1115  sumTZ += scinWeight * zPosition * ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime();
1116 
1117  } // condition of good scin time
1118  ihhit ++;
1119  } // loop over hits of plane
1120  } // loop over planes
1121 
1122  Double_t tmp = sumW * sumZZ - sumZ * sumZ ;
1123  Double_t t0 = ( sumT * sumZZ - sumZ * sumTZ ) / tmp ;
1124  Double_t tmpDenom = sumW * sumTZ - sumZ * sumT;
1125 
1126  if ( TMath::Abs( tmpDenom ) > ( 1 / 10000000000.0 ) ) {
1127 
1128  fBetaNoTrk = tmp / tmpDenom;
1129  fBetaNoTrkChiSq = 0.;
1130  ihhit = 0;
1131 
1132  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){ // Loop over planes
1133  Int_t nphits=fPlanes[ip]->GetNScinHits();
1134  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
1135 
1136  for(Int_t i=0;i<nphits;i++) {
1137  Int_t index=((THcHodoHit*)hodoHits->At(i))->GetPaddleNumber()-1;
1138 
1139  if(twogoodtimes[ihhit]) {
1140 
1141  Double_t zPosition = fPlanes[ip]->GetZpos() + (index%2)*fPlanes[ip]->GetDzpos();
1142  Double_t timeDif = ( ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime() - t0 );
1143 
1144  Double_t sigma = 0.0;
1145  if(fTofUsingInvAdc){
1146  sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoPosSigma[GetScinIndex(ip,index)],2) +
1147  TMath::Power( fHodoNegSigma[GetScinIndex(ip,index)],2) ) );
1148  }
1149  else {
1150  sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoSigmaPos[GetScinIndex(ip,index)],2) +
1151  TMath::Power( fHodoSigmaNeg[GetScinIndex(ip,index)],2) ) );
1152  }
1153 
1154  fBetaNoTrkChiSq += ( ( zPosition / fBetaNoTrk - timeDif ) *
1155  ( zPosition / fBetaNoTrk - timeDif ) ) / ( sigma * sigma );
1156 
1157 
1158  } // condition for good scin time
1159  ihhit++;
1160  } // loop over hits of a plane
1161  } // loop over planes
1162 
1163  Double_t pathNorm = 1.0;
1164 
1165  fBetaNoTrk = fBetaNoTrk * pathNorm;
1166  fBetaNoTrk = fBetaNoTrk / 29.979; // velocity / c
1167 
1168  } // condition for fTmpDenom
1169  else {
1170  fBetaNoTrk = 0.;
1171  fBetaNoTrkChiSq = -2.;
1172  } // else condition for fTmpDenom
1173  //
1175  if ((fNumPlanesBetaCalc==4)&&goodplanetime[0]&&goodplanetime[1]&&goodplanetime[2]&&goodplanetime[3]&&fPlanes[0]->GetNGoodHits()==1&&fPlanes[1]->GetNGoodHits()==1&&fPlanes[2]->GetNGoodHits()==1&&fPlanes[3]->GetNGoodHits()==1) fGoodEventTOFCalib=kTRUE;
1176  if ((fNumPlanesBetaCalc==3)&&goodplanetime[0]&&goodplanetime[1]&&goodplanetime[2]&&fPlanes[0]->GetNGoodHits()==1&&fPlanes[1]->GetNGoodHits()==1&&fPlanes[2]->GetNGoodHits()==1) fGoodEventTOFCalib=kTRUE;
1177  //
1178  //
1179  } else {
1180  fBetaNoTrkChiSq = -10.; // Flag if does not try to find beta
1181  }
1182 }
1183 
1184 //_____________________________________________________________________________
1186 {
1187  return(0);
1188 }
1189 //_____________________________________________________________________________
1191  const ESide side)
1192 {
1193  return(0.0);
1194 }
1195 
1196 
1197 //_____________________________________________________________________________
1199 {
1200 
1201 
1202  Int_t ntracks = tracks.GetLast()+1; // Number of reconstructed tracks
1203  // -------------------------------------------------
1204 
1205  // fDumpOut << " ntrack = " << ntracks << endl;
1206 
1207  if (ntracks > 0 ) {
1208 
1209  // **MAIN LOOP: Loop over all tracks and get corrected time, tof, beta...
1210  vector<Double_t> nPmtHit(ntracks);
1211  vector<Double_t> timeAtFP(ntracks);
1212  fdEdX.reserve(ntracks);
1213  fGoodFlags.reserve(ntracks);
1214  for ( Int_t itrack = 0; itrack < ntracks; itrack++ ) { // Line 133
1215  nPmtHit[itrack]=0;
1216  timeAtFP[itrack]=0;
1217 
1218  THaTrack* theTrack = dynamic_cast<THaTrack*>( tracks.At(itrack) );
1219  if (!theTrack) return -1;
1220 
1221  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
1222  fGoodPlaneTime[ip] = kFALSE;
1223  fNScinHits[ip] = 0;
1224  fNPlaneTime[ip] = 0;
1225  fSumPlaneTime[ip] = 0.;
1226  }
1227  std::vector<Double_t> dedx_temp;
1228  std::vector<std::vector<GoodFlags> > goodflagstmp1;
1229  goodflagstmp1.reserve(fNumPlanesBetaCalc);
1230 #if __cplusplus >= 201103L
1231  fdEdX.push_back(std::move(dedx_temp)); // Create array of dedx per hit
1232  fGoodFlags.push_back(std::move(goodflagstmp1));
1233 #else
1234  fdEdX.push_back(dedx_temp); // Create array of dedx per hit
1235  fGoodFlags.push_back(goodflagstmp1);
1236 #endif
1237  Int_t nFPTime = 0;
1238  Double_t betaChiSq = -3;
1239  Double_t beta = 0;
1240  // timeAtFP[itrack] = 0.;
1241  Double_t sumFPTime = 0.; // Line 138
1242  fNScinHit.push_back(0);
1243 
1244 
1245 
1246  hTime->Reset();
1247  fTOFCalc.clear(); // SAW - Can we
1248  fTOFPInfo.clear(); // SAW - combine these two?
1249  Int_t ihhit = 0; // Hit # overall
1250 
1251  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
1252 
1253  std::vector<GoodFlags> goodflagstmp2;
1254  goodflagstmp2.reserve(fNScinHits[ip]);
1255 #if __cplusplus >= 201103L
1256  fGoodFlags[itrack].push_back(std::move(goodflagstmp2));
1257 #else
1258  fGoodFlags[itrack].push_back(goodflagstmp2);
1259 #endif
1260  fNScinHits[ip] = fPlanes[ip]->GetNScinHits();
1261  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
1262 
1263  Double_t zPos = fPlanes[ip]->GetZpos();
1264  Double_t dzPos = fPlanes[ip]->GetDzpos();
1265 
1266  // first loop over hits with in a single plane
1267  for (Int_t iphit = 0; iphit < fNScinHits[ip]; iphit++ ){
1268  // iphit is hit # within a plane
1269  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(iphit);
1270 
1271  fTOFPInfo.push_back(TOFPInfo());
1272  // Can remove these as we will initialize in the constructor
1273  // fTOFPInfo[ihhit].time_pos = -99.0;
1274  // fTOFPInfo[ihhit].time_neg = -99.0;
1275  // fTOFPInfo[ihhit].keep_pos = kFALSE;
1276  // fTOFPInfo[ihhit].keep_neg = kFALSE;
1277  fTOFPInfo[ihhit].scin_pos_time = 0.0;
1278  fTOFPInfo[ihhit].scin_neg_time = 0.0;
1279  fTOFPInfo[ihhit].hit = hit;
1280  fTOFPInfo[ihhit].planeIndex = ip;
1281  fTOFPInfo[ihhit].hitNumInPlane = iphit;
1282  fTOFPInfo[ihhit].onTrack = kFALSE;
1283 
1284  Int_t paddle = hit->GetPaddleNumber()-1;
1285  Double_t zposition = zPos + (paddle%2)*dzPos;
1286 
1287  Double_t xHitCoord = theTrack->GetX() + theTrack->GetTheta() *
1288  ( zposition ); // Line 183
1289 
1290  Double_t yHitCoord = theTrack->GetY() + theTrack->GetPhi() *
1291  ( zposition ); // Line 184
1292 
1293  Double_t scinTrnsCoord, scinLongCoord;
1294  if ( ( ip == 0 ) || ( ip == 2 ) ){ // !x plane. Line 185
1295  scinTrnsCoord = xHitCoord;
1296  scinLongCoord = yHitCoord;
1297  } else if ( ( ip == 1 ) || ( ip == 3 ) ){ // !y plane. Line 188
1298  scinTrnsCoord = yHitCoord;
1299  scinLongCoord = xHitCoord;
1300  } else { return -1; } // Line 195
1301 
1302  fTOFPInfo[ihhit].scinTrnsCoord = scinTrnsCoord;
1303  fTOFPInfo[ihhit].scinLongCoord = scinLongCoord;
1304 
1305  Double_t scinCenter = fPlanes[ip]->GetPosCenter(paddle) + fPlanes[ip]->GetPosOffset();
1306 
1307  // Index to access the 2d arrays of paddle/scintillator properties
1308  Int_t fPIndex = GetScinIndex(ip,paddle);
1309  Double_t betatrack = theTrack->GetP()/TMath::Sqrt(theTrack->GetP()*theTrack->GetP()+fPartMass*fPartMass);
1310 
1311  if ( TMath::Abs( scinCenter - scinTrnsCoord ) <
1312  ( fPlanes[ip]->GetSize() * 0.5 + fPlanes[ip]->GetHodoSlop() ) ){ // Line 293
1313 
1314  fTOFPInfo[ihhit].onTrack = kTRUE;
1315  Double_t zcor = zposition/(29.979*betatrack)*
1316  TMath::Sqrt(1. + theTrack->GetTheta()*theTrack->GetTheta()
1317  + theTrack->GetPhi()*theTrack->GetPhi());
1318  fTOFPInfo[ihhit].zcor = zcor;
1319  if (fCosmicFlag) {
1320  Double_t zcor = -zposition/(29.979*1.0)*
1321  TMath::Sqrt(1. + theTrack->GetTheta()*theTrack->GetTheta()
1322  + theTrack->GetPhi()*theTrack->GetPhi());
1323  fTOFPInfo[ihhit].zcor = zcor;
1324  }
1325  Double_t tdc_pos = hit->GetPosTDC();
1326  Double_t tdc_neg = hit->GetNegTDC();
1327  //
1328  if( (tdc_pos >=fScinTdcMin && tdc_pos <= fScinTdcMax) &&
1329  (tdc_neg >=fScinTdcMin && tdc_neg <= fScinTdcMax )) {
1330  fTOFPInfo[ihhit].scin_pos_time = hit->GetPosCorrectedTime();
1331  Double_t timep = hit->GetPosCorrectedTime()-zcor;
1332  fTOFPInfo[ihhit].time_pos = timep;
1333  hTime->Fill(timep);
1334  fTOFPInfo[ihhit].scin_neg_time = hit->GetNegCorrectedTime();
1335  Double_t timen = hit->GetNegCorrectedTime()-zcor;
1336  fTOFPInfo[ihhit].time_neg = timen;
1337  hTime->Fill(timen);
1338  } else {
1339  //
1341  if(tdc_pos >=fScinTdcMin && tdc_pos <= fScinTdcMax ) {
1342  Double_t adc_pos = hit->GetPosADC();
1343  Double_t adcamp_pos = hit->GetPosADCpeak();
1344  Double_t pathp = fPlanes[ip]->GetPosLeft() - scinLongCoord;
1345  fTOFPInfo[ihhit].pathp = pathp;
1346  Double_t timep = tdc_pos*fScinTdcToTime;
1347  if(fTofUsingInvAdc) {
1348  timep -= fHodoPosInvAdcOffset[fPIndex]
1349  + pathp/fHodoPosInvAdcLinear[fPIndex]
1350  + fHodoPosInvAdcAdc[fPIndex]
1351  /TMath::Sqrt(TMath::Max(20.0*.020,adc_pos));
1352  } else {
1353  //Double_t tw_corr_pos = fHodoPos_c1[fPIndex]/pow(adcamp_pos/fTdc_Thrs,fHodoPos_c2[fPIndex]) - fHodoPos_c1[fPIndex]/pow(200./fTdc_Thrs, fHodoPos_c2[fPIndex]);
1354  Double_t tw_corr_pos=0.;
1355  pathp=scinLongCoord;
1356  if (adcamp_pos>0) tw_corr_pos = 1./pow(adcamp_pos/fTdc_Thrs,fHodoPos_c2[fPIndex]) - 1./pow(200./fTdc_Thrs, fHodoPos_c2[fPIndex]);
1357  timep += -tw_corr_pos + fHodo_LCoeff[fPIndex]+ pathp/fHodoVelFit[fPIndex];
1358  }
1359  fTOFPInfo[ihhit].scin_pos_time = timep;
1360  timep -= zcor;
1361  fTOFPInfo[ihhit].time_pos = timep;
1362 
1363  hTime->Fill(timep);
1364  }
1365  if(tdc_neg >=fScinTdcMin && tdc_neg <= fScinTdcMax ) {
1366  Double_t adc_neg = hit->GetNegADC();
1367  Double_t adcamp_neg = hit->GetNegADCpeak();
1368  Double_t pathn = scinLongCoord - fPlanes[ip]->GetPosRight();
1369  fTOFPInfo[ihhit].pathn = pathn;
1370  Double_t timen = tdc_neg*fScinTdcToTime;
1371  if(fTofUsingInvAdc) {
1372  timen -= fHodoNegInvAdcOffset[fPIndex]
1373  + pathn/fHodoNegInvAdcLinear[fPIndex]
1374  + fHodoNegInvAdcAdc[fPIndex]
1375  /TMath::Sqrt(TMath::Max(20.0*.020,adc_neg));
1376  } else {
1377  pathn=scinLongCoord ;
1378  Double_t tw_corr_neg =0 ;
1379  if (adcamp_neg >0) tw_corr_neg= 1./pow(adcamp_neg/fTdc_Thrs,fHodoNeg_c2[fPIndex]) - 1./pow(200./fTdc_Thrs, fHodoNeg_c2[fPIndex]);
1380  timen += -tw_corr_neg- 2*fHodoCableFit[fPIndex] + fHodo_LCoeff[fPIndex]- pathn/fHodoVelFit[fPIndex];
1381 
1382  }
1383  fTOFPInfo[ihhit].scin_neg_time = timen;
1384  timen -= zcor;
1385  fTOFPInfo[ihhit].time_neg = timen;
1386  hTime->Fill(timen);
1387  }
1388  } // new fTrackBetaIncludeSinglePmtHits
1389  } // matches else
1390  } // condition for cenetr on a paddle
1391  ihhit++;
1392  } // First loop over hits in a plane <---------
1393 
1394  //-----------------------------------------------------------------------------------------------
1395  //------------- First large loop over scintillator hits ends here --------------------
1396  //-----------------------------------------------------------------------------------------------
1397  }
1398  Int_t nhits=ihhit;
1399 
1400  Double_t TimePeak = DetermineTimePeak(2);
1401  if(TimePeak> 0) {
1402 
1403  for(Int_t ih = 0; ih < nhits; ih++) { // loop over all scintillator hits
1404  if ( ( fTOFPInfo[ih].time_pos > (TimePeak-fTofTolerance) ) && ( fTOFPInfo[ih].time_pos < ( TimePeak + fTofTolerance ) ) ) {
1405  fTOFPInfo[ih].keep_pos=kTRUE;
1406  }
1407  if ( ( fTOFPInfo[ih].time_neg > (TimePeak-fTofTolerance) ) && ( fTOFPInfo[ih].time_neg < ( TimePeak + fTofTolerance ) ) ){
1408  fTOFPInfo[ih].keep_neg=kTRUE;
1409  }
1410  }
1411  }
1412 
1413  //---------------------------------------------------------------------------------------------
1414  // ---------------------- Second loop over scint. hits in a plane -----------------------------
1415  //---------------------------------------------------------------------------------------------
1416 
1417  fdEdX[itrack].reserve(nhits);
1418  fTOFCalc.reserve(nhits);
1419  for(Int_t ih=0; ih < nhits; ih++) {
1420  THcHodoHit *hit = fTOFPInfo[ih].hit;
1421  Int_t iphit = fTOFPInfo[ih].hitNumInPlane;
1422  Int_t ip = fTOFPInfo[ih].planeIndex;
1423  // fDumpOut << " looping over hits = " << ih << " plane = " << ip+1 << endl;
1424  // Flags are used by THcHodoEff
1425  fGoodFlags[itrack][ip].reserve(nhits);
1426  fGoodFlags[itrack][ip].push_back(GoodFlags());
1427  assert( iphit >= 0 && (size_t)iphit < fGoodFlags[itrack][ip].size() );
1428  fGoodFlags[itrack][ip][iphit].onTrack = kFALSE;
1429  fGoodFlags[itrack][ip][iphit].goodScinTime = kFALSE;
1430  fGoodFlags[itrack][ip][iphit].goodTdcNeg = kFALSE;
1431  fGoodFlags[itrack][ip][iphit].goodTdcPos = kFALSE;
1432 
1433  fTOFCalc.push_back(TOFCalc());
1434  // Do we set back to false for each track, or just once per event?
1435  assert( ih >= 0 && (size_t)ih < fTOFCalc.size() );
1436  fTOFCalc[ih].good_scin_time = kFALSE;
1437  // These need a track index too to calculate efficiencies
1438  fTOFCalc[ih].good_tdc_pos = kFALSE;
1439  fTOFCalc[ih].good_tdc_neg = kFALSE;
1440  fTOFCalc[ih].pindex = ip;
1441 
1442  Int_t paddle = hit->GetPaddleNumber()-1;
1443  fTOFCalc[ih].hit_paddle = paddle;
1444  fTOFCalc[ih].good_raw_pad = paddle;
1445 
1446  // Double_t scinCenter = fPlanes[ip]->GetPosCenter(paddle) + fPlanes[ip]->GetPosOffset();
1447  // Double_t scinTrnsCoord = fTOFPInfo[ih].scinTrnsCoord;
1448  // Double_t scinLongCoord = fTOFPInfo[ih].scinLongCoord;
1449 
1450  Int_t fPIndex = GetScinIndex(ip,paddle);
1451 
1452  if (fTOFPInfo[ih].onTrack) {
1453  fGoodFlags[itrack][ip][iphit].onTrack = kTRUE;
1454  if ( fTOFPInfo[ih].keep_pos ) { // 301
1455  fTOFCalc[ih].good_tdc_pos = kTRUE;
1456  fGoodFlags[itrack][ip][iphit].goodTdcPos = kTRUE;
1457  }
1458  if ( fTOFPInfo[ih].keep_neg ) { //
1459  fTOFCalc[ih].good_tdc_neg = kTRUE;
1460  fGoodFlags[itrack][ip][iphit].goodTdcNeg = kTRUE;
1461  }
1462  // ** Calculate ave time for scin and error.
1463  if ( fTOFCalc[ih].good_tdc_pos ){
1464  if ( fTOFCalc[ih].good_tdc_neg ){
1465  fTOFCalc[ih].scin_time = ( fTOFPInfo[ih].scin_pos_time +
1466  fTOFPInfo[ih].scin_neg_time ) / 2.;
1467  fTOFCalc[ih].scin_time_fp = ( fTOFPInfo[ih].time_pos +
1468  fTOFPInfo[ih].time_neg ) / 2.;
1469 
1470  if (fTofUsingInvAdc){
1471  fTOFCalc[ih].scin_sigma = TMath::Sqrt( fHodoPosSigma[fPIndex] * fHodoPosSigma[fPIndex] +
1472  fHodoNegSigma[fPIndex] * fHodoNegSigma[fPIndex] )/2.;
1473  }
1474  else {
1475  fTOFCalc[ih].scin_sigma = TMath::Sqrt( fHodoSigmaPos[fPIndex] * fHodoSigmaPos[fPIndex] +
1476  fHodoSigmaNeg[fPIndex] * fHodoSigmaNeg[fPIndex] )/2.;
1477  }
1478 
1479  fTOFCalc[ih].good_scin_time = kTRUE;
1480  fGoodFlags[itrack][ip][iphit].goodScinTime = kTRUE;
1481  } else{
1482  fTOFCalc[ih].scin_time = fTOFPInfo[ih].scin_pos_time;
1483  fTOFCalc[ih].scin_time_fp = fTOFPInfo[ih].time_pos;
1484 
1485  if (fTofUsingInvAdc){
1486  fTOFCalc[ih].scin_sigma = fHodoPosSigma[fPIndex];
1487  }
1488  else{
1489  fTOFCalc[ih].scin_sigma = fHodoSigmaPos[fPIndex];
1490  }
1491 
1492  fTOFCalc[ih].good_scin_time = kTRUE;
1493  fGoodFlags[itrack][ip][iphit].goodScinTime = kTRUE;
1494  }
1495  } else {
1496  if ( fTOFCalc[ih].good_tdc_neg ){
1497  fTOFCalc[ih].scin_time = fTOFPInfo[ih].scin_neg_time;
1498  fTOFCalc[ih].scin_time_fp = fTOFPInfo[ih].time_neg;
1499  if (fTofUsingInvAdc){
1500  fTOFCalc[ih].scin_sigma = fHodoNegSigma[fPIndex];
1501  }
1502  else{
1503  fTOFCalc[ih].scin_sigma = fHodoSigmaNeg[fPIndex];
1504  }
1505  fTOFCalc[ih].good_scin_time = kTRUE;
1506  fGoodFlags[itrack][ip][iphit].goodScinTime = kTRUE;
1507  }
1508  } // In h_tof.f this includes the following if condition for time at focal plane
1509  // // because it is written in FORTRAN code
1510 
1511  // c Get time at focal plane
1512  if ( fTOFCalc[ih].good_scin_time ){
1513 
1514  // scin_time_fp doesn't need to be an array
1515  // Is this any different than the average of time_pos and time_neg?
1516  // Double_t scin_time_fp = ( fTOFPInfo[ih].time_pos +
1517  // fTOFPInfo[ih].time_neg ) / 2.;
1518  Double_t scin_time_fp = fTOFCalc[ih].scin_time_fp;
1519 
1520  sumFPTime = sumFPTime + scin_time_fp;
1521  nFPTime ++;
1522 
1523  fSumPlaneTime[ip] = fSumPlaneTime[ip] + scin_time_fp;
1524  fNPlaneTime[ip] ++;
1525  fNScinHit[itrack] ++;
1526 
1527  if ( ( fTOFCalc[ih].good_tdc_pos ) && ( fTOFCalc[ih].good_tdc_neg ) ){
1528  nPmtHit[itrack] = nPmtHit[itrack] + 2;
1529  } else {
1530  nPmtHit[itrack] = nPmtHit[itrack] + 1;
1531  }
1532 
1533  fdEdX[itrack].push_back(0.0);
1534  assert( fNScinHit[itrack] > 0 && (size_t)fNScinHit[itrack] < fdEdX[itrack].size()+1 );
1535 
1536  // --------------------------------------------------------------------------------------------
1537  if ( fTOFCalc[ih].good_tdc_pos ){
1538  if ( fTOFCalc[ih].good_tdc_neg ){
1539  fdEdX[itrack][fNScinHit[itrack]-1]=
1540  TMath::Sqrt( TMath::Max( 0., hit->GetPosADC() * hit->GetNegADC() ) );
1541  } else{
1542  fdEdX[itrack][fNScinHit[itrack]-1]=
1543  TMath::Max( 0., hit->GetPosADC() );
1544  }
1545  } else{
1546  if ( fTOFCalc[ih].good_tdc_neg ){
1547  fdEdX[itrack][fNScinHit[itrack]-1]=
1548  TMath::Max( 0., hit->GetNegADC() );
1549  } else{
1550  fdEdX[itrack][fNScinHit[itrack]-1]=0.0;
1551  }
1552  }
1553  // --------------------------------------------------------------------------------------------
1554 
1555 
1556  } // time at focal plane condition
1557  } // on track condition
1558 
1559  // ** See if there are any good time measurements in the plane.
1560  if ( fTOFCalc[ih].good_scin_time ){
1561  fGoodPlaneTime[ip] = kTRUE;
1562  fTOFCalc[ih].dedx = fdEdX[itrack][fNScinHit[itrack]-1];
1563  } else {
1564  fTOFCalc[ih].dedx = 0.0;
1565  }
1566 
1567  } // Second loop over hits of a scintillator plane ends here
1568  theTrack->SetGoodPlane3( fGoodPlaneTime[2] ? 1 : 0 );
1569  if (fNumPlanesBetaCalc==4) theTrack->SetGoodPlane4( fGoodPlaneTime[3] ? 1 : 0 );
1570  //
1571  //------------------------------------------------------------------------------
1572  //------------------------------------------------------------------------------
1573  //------------------------------------------------------------------------------
1574  //------------------------------------------------------------------------------
1575  //------------------------------------------------------------------------------
1576  //------------------------------------------------------------------------------
1577  //------------------------------------------------------------------------------
1578  //------------------------------------------------------------------------------
1579 
1580  // * * Fit beta if there are enough time measurements (one upper, one lower)
1581  // From h_tof_fit
1582  if ( ( ( fGoodPlaneTime[0] ) || ( fGoodPlaneTime[1] ) ) &&
1583  ( ( fGoodPlaneTime[2] ) || ( fGoodPlaneTime[3] ) ) ){
1584 
1585  Double_t sumW = 0.;
1586  Double_t sumT = 0.;
1587  Double_t sumZ = 0.;
1588  Double_t sumZZ = 0.;
1589  Double_t sumTZ = 0.;
1590 
1591  for(Int_t ih=0; ih < nhits; ih++) {
1592  Int_t ip = fTOFPInfo[ih].planeIndex;
1593 
1594  if ( fTOFCalc[ih].good_scin_time ) {
1595 
1596  Double_t scinWeight = 1 / ( fTOFCalc[ih].scin_sigma * fTOFCalc[ih].scin_sigma );
1597  Double_t zPosition = ( fPlanes[ip]->GetZpos()
1598  +( fTOFCalc[ih].hit_paddle % 2 ) *
1599  fPlanes[ip]->GetDzpos() );
1600 
1601  sumW += scinWeight;
1602  sumT += scinWeight * fTOFCalc[ih].scin_time;
1603  sumZ += scinWeight * zPosition;
1604  sumZZ += scinWeight * ( zPosition * zPosition );
1605  sumTZ += scinWeight * zPosition * fTOFCalc[ih].scin_time;
1606 
1607  } // condition of good scin time
1608  } // loop over hits
1609 
1610  Double_t tmp = sumW * sumZZ - sumZ * sumZ ;
1611  Double_t t0 = ( sumT * sumZZ - sumZ * sumTZ ) / tmp ;
1612  Double_t tmpDenom = sumW * sumTZ - sumZ * sumT;
1613 
1614  if ( TMath::Abs( tmpDenom ) > ( 1 / 10000000000.0 ) ) {
1615 
1616  beta = tmp / tmpDenom;
1617  betaChiSq = 0.;
1618 
1619  for(Int_t ih=0; ih < nhits; ih++) {
1620  Int_t ip = fTOFPInfo[ih].planeIndex;
1621 
1622  if ( fTOFCalc[ih].good_scin_time ){
1623 
1624  Double_t zPosition = ( fPlanes[ip]->GetZpos() + ( fTOFCalc[ih].hit_paddle % 2 ) *
1625  fPlanes[ip]->GetDzpos() );
1626  Double_t timeDif = ( fTOFCalc[ih].scin_time - t0 );
1627  betaChiSq += ( ( zPosition / beta - timeDif ) *
1628  ( zPosition / beta - timeDif ) ) /
1629  ( fTOFCalc[ih].scin_sigma * fTOFCalc[ih].scin_sigma );
1630 
1631  } // condition for good scin time
1632  } // loop over hits
1633 
1634  Double_t pathNorm = TMath::Sqrt( 1. + theTrack->GetTheta() * theTrack->GetTheta() +
1635  theTrack->GetPhi() * theTrack->GetPhi() );
1636  // Take angle into account
1637  beta = beta / pathNorm;
1638  beta = beta / 29.979; // velocity / c
1639 
1640  } // condition for fTmpDenom
1641  else {
1642  beta = 0.;
1643  betaChiSq = -2.;
1644  } // else condition for fTmpDenom
1645  } else {
1646  beta = 0.;
1647  betaChiSq = -1;
1648  }
1649 
1650  if ( nFPTime != 0 ){
1651  timeAtFP[itrack] = ( sumFPTime / nFPTime );
1652  }
1653  //
1654  // ---------------------------------------------------------------------------
1655 
1656  Double_t FPTimeSum=0.0;
1657  Int_t nFPTimeSum=0;
1658  Int_t nGoodPlanesHit=0;
1659  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
1660  if ( fNPlaneTime[ip] != 0 ){
1661  nGoodPlanesHit++;
1662  fFPTime[ip] = ( fSumPlaneTime[ip] / fNPlaneTime[ip] );
1663  FPTimeSum += fSumPlaneTime[ip];
1664  nFPTimeSum += fNPlaneTime[ip];
1665  } else {
1666  fFPTime[ip] = 1000. * ( ip + 1 );
1667  }
1668  }
1669  Double_t fptime=-2000;
1670  fptime=fStartTime;
1671  if (nGoodPlanesHit>=3) fptime = FPTimeSum/nFPTimeSum;
1672  fFPTimeAll = fptime;
1673  Double_t dedx=0.0;
1674  for(UInt_t ih=0;ih<fTOFCalc.size();ih++) {
1675  if(fTOFCalc[ih].good_scin_time) {
1676  dedx = fTOFCalc[ih].dedx;
1677  break;
1678  }
1679  }
1680  theTrack->SetDedx(dedx);
1681  theTrack->SetFPTime(fptime);
1682  theTrack->SetBeta(beta);
1683  theTrack->SetBetaChi2( betaChiSq );
1684  theTrack->SetNPMT(nPmtHit[itrack]);
1685 
1686 
1687  } // Main loop over tracks ends here.
1688 
1689  } // If condition for at least one track
1690 
1691 
1692  //OriginalTrackEffTest();
1693  TrackEffTest();
1694  //
1695  CalcCluster();
1696 
1697  return 0;
1698 
1699 }
1700 //
1702 {
1703  // THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
1704  // cout << app->GetName() << endl;
1705  const Int_t MaxNCluster=5;
1706  std::vector<Int_t > iw(MaxNCluster,0);
1707  std::vector<Double_t > dw(MaxNCluster,0);
1708  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
1709  fNCluster.push_back(0);
1710  fClusterSize.push_back(iw);
1711  fClusterXPos.push_back(dw);
1712  fClusterYPos.push_back(dw);
1713  }
1714  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
1715  Double_t pl_xypos=0;
1716  Double_t pl_calcpos=0;
1717  Double_t pl_zpos=0;
1718  Int_t num_good_pad=0;
1719  Double_t pl_x=0,pl_y=0;
1720  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
1721  Int_t prev_padnum=-100;
1722  for (Int_t iphit = 0; iphit < fPlanes[ip]->GetNScinHits(); iphit++ ){
1723  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(iphit);
1724  if ( hit->GetTwoGoodTimes() ) {
1725  Int_t padind = hit->GetPaddleNumber()-1;
1726  Int_t padnum = padind+1;
1727  if (ip==0 || ip==2) pl_x = fPlanes[ip]->GetPosCenter(padind)+ fPlanes[ip]->GetPosOffset();
1728  if (ip==0 || ip==2) pl_y = ((THcHodoHit*)hodoHits->At(iphit))->GetCalcPosition();
1729  if (ip==1 || ip==3) pl_y = fPlanes[ip]->GetPosCenter(padind)+ fPlanes[ip]->GetPosOffset();
1730  if (ip==1 || ip==3) pl_x = ((THcHodoHit*)hodoHits->At(iphit))->GetCalcPosition();
1731  pl_xypos+=fPlanes[ip]->GetPosCenter(padind)+ fPlanes[ip]->GetPosOffset();
1732  pl_calcpos=((THcHodoHit*)hodoHits->At(iphit))->GetCalcPosition();
1733  pl_zpos+=fPlanes[ip]->GetZpos()+ (padind%2)*fPlanes[ip]->GetDzpos();
1734  num_good_pad++;
1735  if ( fNCluster[ip]>0 && abs(padnum-prev_padnum)==1 && fClusterSize[ip][fNCluster[ip]-1]==1) {
1736  fClusterSize[ip][fNCluster[ip]-1]=fClusterSize[ip][fNCluster[ip]-1]+1;
1737  fClusterXPos[ip][fNCluster[ip]-1]+=pl_x;
1738  fClusterYPos[ip][fNCluster[ip]-1]+=pl_y;
1739  // cout << "Add to cluster pl = " << ip+1 << " hit = " << iphit << " pad = " << padnum << " clus = " << fNCluster[ip] << " cl size = " << fClusterSize[ip][fNCluster[ip]-1] << " Xpos " << pl_x << " Ypos = " << pl_y << " postime = " << hit->GetPosTOFCorrectedTime() << " negtime = " << hit->GetNegTOFCorrectedTime() << endl;
1740  } else {
1741  if (fNCluster[ip]<MaxNCluster) fNCluster[ip]++;
1742  fClusterSize[ip][fNCluster[ip]-1]=1;
1743  fClusterXPos[ip][fNCluster[ip]-1]=pl_x;
1744  fClusterYPos[ip][fNCluster[ip]-1]=pl_y;
1745  // cout << " New clus pl = " << ip+1 << " hit = " << iphit << " pad = " << padnum << " clus = " << fNCluster[ip] << " cl size = " << fClusterSize[ip][fNCluster[ip]-1] << " Xpos = " << pl_x << " Ypos = " << pl_y << " postime = " << hit->GetPosTOFCorrectedTime() << " negtime = " << hit->GetNegTOFCorrectedTime() << endl;
1746  }
1747  prev_padnum=padnum;
1748  }
1749  }
1750  //
1751  for (Int_t ic = 0; ic < fNCluster[ip]; ic++ ){
1752  fClusterXPos[ip][ic]/=fClusterSize[ip][ic];
1753  fClusterYPos[ip][ic]/=fClusterSize[ip][ic];
1754  // cout << " Cluster = " << ic+1 << " Xpos = " << fClusterXPos[ip][ic] << " Ypos = " << fClusterYPos[ip][ic] << endl;
1755  }
1756  //
1757  if (num_good_pad !=0 ) {
1758  pl_xypos=pl_xypos/num_good_pad;
1759  pl_calcpos=pl_calcpos/num_good_pad;
1760  pl_zpos=pl_zpos/num_good_pad;
1761  } else {
1762  pl_xypos = kBig;
1763  pl_calcpos = kBig;
1764  pl_zpos = kBig;
1765  }
1766  // fPlanes[ip]->SetTranversePos(pl_xypos);
1767  //fPlanes[ip]->SetLongPos(pl_calcpos);
1768  }
1769  //
1770  // analyze clusters
1771  Int_t best_cluster[4]={-1,-1,-1,-1};
1772  Double_t diffx_test;
1773  Double_t diffx;
1774  Double_t diffy_test;
1775  Double_t diffy;
1776  Int_t pl1,pl2;
1777  for (Int_t nch = 0; nch < 2; nch++ ){
1778  if (nch==0) pl1=0;
1779  if (nch==1) pl1=2;
1780  pl2=pl1+1;
1781  diffx_test=1000;
1782  diffy_test=1000;
1783  if ( fNCluster[pl1]>=1 && fNCluster[pl2]>=1 ) {
1784  for (Int_t ic1 = 0; ic1 < fNCluster[pl1]; ic1++ ){
1785  for (Int_t ic2 = 0; ic2 < fNCluster[pl2]; ic2++ ){
1786  diffx= abs(fClusterXPos[pl1][ic1]-fClusterXPos[pl2][ic2]);
1787  diffy= abs(fClusterYPos[pl1][ic1]-fClusterYPos[pl2][ic2]);
1788  if ( (ic1==0 && ic2==0) || (diffx <=diffx_test && diffy <=diffy_test)) {
1789  diffx_test=diffx;
1790  diffy_test=diffy;
1791  best_cluster[pl1]=ic1;
1792  best_cluster[pl2]=ic2;
1793  }
1794  }}
1795  } else {
1796  if (fNCluster[pl1]==1) best_cluster[pl1]=0;
1797  if (fNCluster[pl2]==1) best_cluster[pl2]=0;
1798  }
1799  }
1800  //
1801  Int_t pl_test1[4]={0,1,2,3};
1802  Int_t pl_test2[4]={2,3,0,1};
1803  for (Int_t npl = 0; npl < 4; npl++ ){
1804  pl1=pl_test1[npl];
1805  pl2=pl_test2[npl];
1806  if (fNCluster[pl1]>0 && best_cluster[pl1]==-1 && fNCluster[pl2]>0 && best_cluster[pl2]>-1) {
1807  diffx_test=1000;
1808  diffy_test=1000;
1809  for (Int_t ic1 = 0; ic1 < fNCluster[pl1]; ic1++ ){
1810  diffx= abs(fClusterXPos[pl1][ic1]-fClusterXPos[pl2][best_cluster[pl2]]);
1811  diffy= abs(fClusterYPos[pl1][ic1]-fClusterYPos[pl2][best_cluster[pl2]]);
1812  if ( (diffx <=diffx_test && diffy <=diffy_test)) {
1813  diffx_test=diffx;
1814  diffy_test=diffy;
1815  best_cluster[pl1]=ic1;
1816  }
1817  }
1818  }
1819  }
1820  //
1821  //
1822 
1823  //
1824  for (Int_t npl = 0; npl < 4; npl++ ){
1825  /*
1826  if (best_cluster[npl]==-1) {
1827  cout << " PLane = " << npl+1 << " no best cluster " << endl;
1828  } else {
1829  cout << " plane = " << npl+1 << " xpos = " << fClusterXPos[npl][best_cluster[npl]] << " ypos = " << fClusterYPos[npl][best_cluster[npl]] << endl;
1830  }
1831  */
1832  if (best_cluster[npl]!=-1) fPlanes[npl]->SetScinYPos( fClusterYPos[npl][best_cluster[npl]] );
1833  if (best_cluster[npl]!=-1) fPlanes[npl]->SetScinXPos( fClusterXPos[npl][best_cluster[npl]] );
1834  }
1835  //
1836 }
1837 //
1839 {
1840  Double_t PadLow[4];
1841  Double_t PadHigh[4];
1842  // assume X planes are 0,2 and Y planes are 1,3
1843  PadLow[0]=fxLoScin[0];
1844  PadLow[2]=fxLoScin[1];
1845  PadLow[1]=fyLoScin[0];
1846  PadLow[3]=fyLoScin[1];
1847  PadHigh[0]=fxHiScin[0];
1848  PadHigh[2]=fxHiScin[1];
1849  PadHigh[1]=fyHiScin[0];
1850  PadHigh[3]=fyHiScin[1];
1851  //
1852  Bool_t efftest_debug = kFALSE;
1853  if (efftest_debug) cout << " spec = " << GetApparatus()->GetName()[0] << endl;
1854  Double_t PadPosLo[4];
1855  Double_t PadPosHi[4];
1856  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
1857  Double_t lowtemp=fPlanes[ip]->GetPosCenter(PadLow[ip]-1)+ fPlanes[ip]->GetPosOffset();
1858  Double_t hitemp=fPlanes[ip]->GetPosCenter(PadHigh[ip]-1)+ fPlanes[ip]->GetPosOffset();
1859  if (lowtemp < hitemp) {
1860  PadPosLo[ip]=lowtemp;
1861  PadPosHi[ip]=hitemp;
1862  } else {
1863  PadPosLo[ip]=hitemp;
1864  PadPosHi[ip]=lowtemp;
1865  }
1866  }
1867  //
1868  const Int_t MaxNClus=5;
1869  std::vector<Int_t > iw(MaxNClus,0);
1870  std::vector<Double_t > dw(MaxNClus,0);
1871  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
1872  fNClust.push_back(0);
1873  fClustSize.push_back(iw);
1874  fClustPos.push_back(dw);
1875  }
1876  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
1877  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
1878  Int_t prev_padnum=-100;
1879  for (Int_t iphit = 0; iphit < fPlanes[ip]->GetNScinHits(); iphit++ ){
1880  THcHodoHit *hit = (THcHodoHit*)hodoHits->At(iphit);
1881  Int_t padnum = hit->GetPaddleNumber();
1882  if ( hit->GetTwoGoodTimes() ) {
1883  if ( padnum==prev_padnum+1 ) {
1884  fClustSize[ip][fNClust[ip]-1]=fClustSize[ip][fNClust[ip]-1]+1;
1885  fClustPos[ip][fNClust[ip]-1]+=fPlanes[ip]->GetPosCenter(padnum-1)+ fPlanes[ip]->GetPosOffset();
1886  if (efftest_debug) cout << "Add to cluster pl = " << ip+1 << " hit = " << iphit << " pad = " << padnum << " clus = " << fNClust[ip] << " cl size = " << fClustSize[ip][fNClust[ip]-1] << " pos " << fPlanes[ip]->GetPosCenter(padnum-1)+ fPlanes[ip]->GetPosOffset() << endl;
1887  } else {
1888  if (fNClust[ip]<MaxNClus) fNClust[ip]++;
1889  fClustSize[ip][fNClust[ip]-1]=1;
1890  fClustPos[ip][fNClust[ip]-1]=fPlanes[ip]->GetPosCenter(padnum-1)+ fPlanes[ip]->GetPosOffset();
1891  if (efftest_debug) cout << " New clus pl = " << ip+1 << " hit = " << iphit << " pad = " << padnum << " clus = " << fNClust[ip] << " cl size = " << fClustSize[ip][fNClust[ip]-1] << " pos " << fPlanes[ip]->GetPosCenter(padnum-1)+ fPlanes[ip]->GetPosOffset() << endl;
1892  }
1893  prev_padnum=padnum;
1894  }
1895  if (!(hit->GetTwoGoodTimes()) && efftest_debug) cout << "no two good times plane = " << ip+1 << " hit = " << iphit << endl;
1896  }
1897  }
1898  //
1899  Bool_t inside_bound[4][MaxNClus];
1900  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
1901  fPlanes[ip]->SetNumberClusters(fNClust[ip]);
1902  for(Int_t ic = 0; ic <fNClust[ip] ; ic++ ) {
1903  fClustPos[ip][ic]=fClustPos[ip][ic]/fClustSize[ip][ic];
1904  fPlanes[ip]->SetCluster(ic,fClustPos[ip][ic]);
1905  fPlanes[ip]->SetClusterSize(ic,fClustSize[ip][ic]);
1906  inside_bound[ip][ic] = fClustPos[ip][ic]>=PadPosLo[ip] && fClustPos[ip][ic]<=PadPosHi[ip];
1907  if (efftest_debug) cout << "plane = " << ip+1 << " Cluster = " << ic+1 << " size = " << fClustSize[ip][ic]<< " pos = " << fClustPos[ip][ic] << " inside = " << inside_bound[ip][ic] << " lo = " << PadPosLo[ip]<< " hi = " << PadPosHi[ip]<< endl;
1908  }
1909  }
1910  //
1911  Int_t MaxClusterSize=3;
1912  Int_t good_for_track_test[4][MaxNClus];
1913  Int_t sum_good_track_test[4]={0,0,0,0};
1914  Int_t num_good_plane_hit=0;
1915  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
1916  for(Int_t ic = 0; ic <fNClust[ip] ; ic++ ) {
1917  if (inside_bound[ip][ic] && fClustSize[ip][ic]<=MaxClusterSize) {
1918  fPlanes[ip]->SetClusterFlag(ic,1.);
1919  good_for_track_test[ip][ic]=1;
1920  sum_good_track_test[ip]++;
1921  if (sum_good_track_test[ip]==1) num_good_plane_hit++;
1922  } else {
1923  good_for_track_test[ip][ic]=0;
1924  }
1925  if (efftest_debug) cout << " ip " << ip+1 << " clus = " << ic << " good for track = " << good_for_track_test[ip][ic] << endl;
1926  }
1927  if (efftest_debug) cout << " ip = " << ip+1 << " sum_good_track_test = " << sum_good_track_test[ip] << endl;
1928  }
1929  if (efftest_debug) cout << " number of planes hits = " << num_good_plane_hit << endl;
1930  //
1931  Bool_t xdiffTest=kFALSE;
1932  Bool_t ydiffTest=kFALSE;
1933  fGoodScinHits = 0;
1934  if (efftest_debug) cout << " fTrackEffTestNScinPlanes = " << fTrackEffTestNScinPlanes << endl;
1935  if ( (fTrackEffTestNScinPlanes == 4 || fTrackEffTestNScinPlanes == 3) && num_good_plane_hit==4) {
1936 
1937  // check for matching clusters in the X planes assumed to be planes 0 and 2
1938  for(Int_t ic0 = 0; ic0 <fNClust[0] ; ic0++ ) {
1939  for(Int_t ic2 = 0; ic2 <fNClust[2] ; ic2++ ) {
1940  if (good_for_track_test[0][ic0] && good_for_track_test[2][ic2]) {
1941  Double_t x1_proj = fClustPos[0][ic0]*(1+fRatio_xpfp_to_xfp*(fPlanes[2]->GetZpos()-fPlanes[0]->GetZpos())); // project X1 to X2 Z position
1942  xdiffTest= TMath::Abs(x1_proj-fClustPos[2][ic2])<trackeff_scint_xdiff_max;
1943  if (xdiffTest) fPlanes[0]->SetClusterUsedFlag(ic0,1.);
1944  if (xdiffTest) fPlanes[2]->SetClusterUsedFlag(ic2,1.);
1945  }
1946  }
1947  }
1948  // check for matching clusters in the Y planes assumed to be planes 1 and 3
1949  for(Int_t ic1 = 0; ic1 <fNClust[1] ; ic1++ ) {
1950  for(Int_t ic3 = 0; ic3 <fNClust[3] ; ic3++ ) {
1951  if (good_for_track_test[1][ic1] && good_for_track_test[3][ic3]) {
1952  ydiffTest= TMath::Abs(fClustPos[1][ic1]-fClustPos[3][ic3])<trackeff_scint_ydiff_max;
1953  if (ydiffTest) fPlanes[1]->SetClusterUsedFlag(ic1,1.);
1954  if (ydiffTest) fPlanes[3]->SetClusterUsedFlag(ic3,1.);
1955  }
1956  }
1957  }
1958  if (xdiffTest && ydiffTest) fGoodScinHits = 1;
1959  if (efftest_debug) cout << " 4 good planes xdiff = " << xdiffTest << " ydiff = " << ydiffTest << endl;
1960  }
1961  //
1962  if (fTrackEffTestNScinPlanes == 3 && num_good_plane_hit==3) {
1963  xdiffTest=kFALSE;
1964  ydiffTest=kFALSE;
1965  // Check if two X planes hit
1966  if (sum_good_track_test[0]>0&&sum_good_track_test[2]>0) {
1967  for(Int_t ic0 = 0; ic0 <fNClust[0] ; ic0++ ) {
1968  for(Int_t ic2 = 0; ic2 <fNClust[2] ; ic2++ ) {
1969  if (good_for_track_test[0][ic0] && good_for_track_test[2][ic2]) {
1970  xdiffTest= TMath::Abs(fClustPos[0][ic0]-fClustPos[2][ic2])<trackeff_scint_xdiff_max;
1971  }
1972  }
1973  }
1974  ydiffTest = kTRUE;
1975  }
1976  // Check if two Y planes hit
1977  if ((sum_good_track_test[1]>0||sum_good_track_test[3]>0)) {
1978  for(Int_t ic1 = 0; ic1 <fNClust[1] ; ic1++ ) {
1979  for(Int_t ic3 = 0; ic3 <fNClust[3] ; ic3++ ) {
1980  if (good_for_track_test[1][ic1] && good_for_track_test[3][ic3]) {
1981  ydiffTest= TMath::Abs(fClustPos[1][ic1]-fClustPos[3][ic3])<trackeff_scint_ydiff_max;
1982  }
1983  }
1984  xdiffTest = kTRUE;
1985  }
1986  }
1987  if (xdiffTest && ydiffTest) fGoodScinHits = 1;
1988  if (efftest_debug) cout << " 3 good planes xdiff = " << xdiffTest << " ydiff = " << ydiffTest << endl;
1989  }
1990  if (efftest_debug) cout << " ************" << endl;
1991  //
1992 }
1993 //
1994 //
1995 //
1997 {
2002  //************************now look at some hodoscope tests
2003  // *second, we move the scintillators. here we use scintillator cuts to see
2004  // *if a track should have been found.
2005  cout << " enter track eff" << fNumPlanesBetaCalc << endl;
2006  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
2007  cout << " loop over planes " << ip+1 << endl;
2008 
2009  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
2010  // TClonesArray* scinPosTDC = fPlanes[ip]->GetPosTDC();
2011  // TClonesArray* scinNegTDC = fPlanes[ip]->GetNegTDC();
2012 
2013  fNScinHits[ip] = fPlanes[ip]->GetNScinHits();
2014  cout << " hits = " << fNScinHits[ip] << endl;
2015  for (Int_t iphit = 0; iphit < fNScinHits[ip]; iphit++ ){
2016  Int_t paddle = ((THcHodoHit*)hodoHits->At(iphit))->GetPaddleNumber()-1;
2017 
2018  fScinHitPaddle[ip][paddle] = 1;
2019  cout << " hit = " << iphit+1 << " " << paddle+1 << endl;
2020  }
2021  }
2022 
2023  // *next, look for clusters of hits in a scin plane. a cluster is a group of
2024  // *adjacent scintillator hits separated by a non-firing scintillator.
2025  // *Wwe count the number of three adjacent scintillators too. (A signle track
2026  // *shouldn't fire three adjacent scintillators.
2027 
2028  for(Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ) {
2029  // Planes ip = 0 = 1X
2030  // Planes ip = 2 = 2X
2031  fNClust.push_back(0);
2032  fThreeScin.push_back(0);
2033  }
2034 
2035  // *look for clusters in x planes... (16 scins) !this assume both x planes have same
2036  // *number of scintillators.
2037  cout << " looking for cluster in x planes" << endl;
2038  Int_t icount;
2039  for (Int_t ip = 0; ip < 3; ip +=2 ) {
2040  icount = 0;
2041  if ( fScinHitPaddle[ip][0] == 1 ) icount ++;
2042  cout << "plane =" << ip << "check if paddle 1 hit " << icount << endl;
2043 
2044  for (Int_t ipaddle = 0; ipaddle < (Int_t) fNPaddle[ip] - 1; ipaddle++ ){
2045  // !look for number of clusters of 1 or more hits
2046  if ( ( fScinHitPaddle[ip][ipaddle] == 0 ) &&
2047  ( fScinHitPaddle[ip][ipaddle + 1] == 1 ) )
2048  icount ++;
2049  cout << " paddle = " << ipaddle+1 << " " << icount << endl;
2050 
2051  } // Loop over paddles
2052  cout << "Two cluster in plane = " << ip+1 << " " << icount << endl;
2053  fNClust[ip] = icount;
2054  icount = 0;
2055 
2056  for (Int_t ipaddle = 0; ipaddle < (Int_t) fNPaddle[ip] - 2; ipaddle++ ){
2057  // !look for three or more adjacent hits
2058 
2059  if ( ( fScinHitPaddle[ip][ipaddle] == 1 ) &&
2060  ( fScinHitPaddle[ip][ipaddle + 1] == 1 ) &&
2061  ( fScinHitPaddle[ip][ipaddle + 2] == 1 ) )
2062  icount ++;
2063  } // Second loop over paddles
2064  cout << "Three clusters in plane = " << ip+1 << " " << icount << endl;
2065 
2066  if ( icount > 0 )
2067  fThreeScin[ip] = 1;
2068 
2069  } // Loop over X plane
2070 
2071  // *look for clusters in y planes... (10 scins) !this assume both y planes have same
2072  // *number of scintillators.
2073  cout << " looking for cluster in y planes" << endl;
2074 
2075  for (Int_t ip = 1; ip < fNumPlanesBetaCalc; ip +=2 ) {
2076  // Planes ip = 1 = 1Y
2077  // Planes ip = 3 = 2Y
2078 
2079  icount = 0;
2080  if ( fScinHitPaddle[ip][0] == 1 ) icount ++;
2081  cout << "plane =" << ip << "check if paddle 1 hit " << icount << endl;
2082 
2083  for (Int_t ipaddle = 0; ipaddle < (Int_t) fNPaddle[ip] - 1; ipaddle++ ){
2084  // !look for number of clusters of 1 or more hits
2085 
2086  if ( ( fScinHitPaddle[ip][ipaddle] == 0 ) &&
2087  ( fScinHitPaddle[ip][ipaddle + 1] == 1 ) )
2088  icount ++;
2089  cout << " paddle = " << ipaddle+1 << " " << icount << endl;
2090 
2091  } // Loop over Y paddles
2092  cout << "Two cluster in plane = " << ip+1 << " " << icount << endl;
2093 
2094  fNClust[ip] = icount;
2095  icount = 0;
2096 
2097  for (Int_t ipaddle = 0; ipaddle < (Int_t) fNPaddle[ip] - 2; ipaddle++ ){
2098  // !look for three or more adjacent hits
2099 
2100  if ( ( fScinHitPaddle[ip][ipaddle] == 1 ) &&
2101  ( fScinHitPaddle[ip][ipaddle + 1] == 1 ) &&
2102  ( fScinHitPaddle[ip][ipaddle + 2] == 1 ) )
2103  icount ++;
2104 
2105  } // Second loop over Y paddles
2106  cout << "Three clusters in plane = " << ip+1 << " " << icount << endl;
2107 
2108  if ( icount > 0 )
2109  fThreeScin[ip] = 1;
2110 
2111  }// Loop over Y planes
2112 
2113 
2114  // *next we mask out the edge scintillators, and look at triggers that happened
2115  // *at the center of the acceptance. To change which scins are in the mask
2116  // *change the values of h*loscin and h*hiscin in htracking.param
2117 
2118  // fGoodScinHits = 0;
2119  for (Int_t ifidx = fxLoScin[0]; ifidx < (Int_t) fxHiScin[0]; ifidx ++ ){
2120  fGoodScinHitsX.push_back(0);
2121  }
2122 
2123  fHitSweet1X=0;
2124  fHitSweet2X=0;
2125  fHitSweet1Y=0;
2126  fHitSweet2Y=0;
2127  // *first x plane. first see if there are hits inside the scin region
2128  for (Int_t ifidx = fxLoScin[0]-1; ifidx < fxHiScin[0]; ifidx ++ ){
2129  if ( fScinHitPaddle[0][ifidx] == 1 ){
2130  fHitSweet1X = 1;
2131  fSweet1XScin = ifidx + 1;
2132  }
2133  }
2134 
2135  // * next make sure nothing fired outside the good region
2136  for (Int_t ifidx = 0; ifidx < fxLoScin[0]-1; ifidx ++ ){
2137  if ( fScinHitPaddle[0][ifidx] == 1 ){ fHitSweet1X = -1; }
2138  }
2139  for (Int_t ifidx = fxHiScin[0]; ifidx < (Int_t) fNPaddle[0]; ifidx ++ ){
2140  if ( fScinHitPaddle[0][ifidx] == 1 ){ fHitSweet1X = -1; }
2141  }
2142 
2143  // *second x plane. first see if there are hits inside the scin region
2144  for (Int_t ifidx = fxLoScin[1]-1; ifidx < fxHiScin[1]; ifidx ++ ){
2145  if ( fScinHitPaddle[2][ifidx] == 1 ){
2146  fHitSweet2X = 1;
2147  fSweet2XScin = ifidx + 1;
2148  }
2149  }
2150  // * next make sure nothing fired outside the good region
2151  for (Int_t ifidx = 0; ifidx < fxLoScin[1]-1; ifidx ++ ){
2152  if ( fScinHitPaddle[2][ifidx] == 1 ){ fHitSweet2X = -1; }
2153  }
2154  for (Int_t ifidx = fxHiScin[1]; ifidx < (Int_t) fNPaddle[2]; ifidx ++ ){
2155  if ( fScinHitPaddle[2][ifidx] == 1 ){ fHitSweet2X = -1; }
2156  }
2157 
2158  // *first y plane. first see if there are hits inside the scin region
2159  for (Int_t ifidx = fyLoScin[0]-1; ifidx < fyHiScin[0]; ifidx ++ ){
2160  if ( fScinHitPaddle[1][ifidx] == 1 ){
2161  fHitSweet1Y = 1;
2162  fSweet1YScin = ifidx + 1;
2163  }
2164  }
2165  // * next make sure nothing fired outside the good region
2166  for (Int_t ifidx = 0; ifidx < fyLoScin[0]-1; ifidx ++ ){
2167  if ( fScinHitPaddle[1][ifidx] == 1 ){ fHitSweet1Y = -1; }
2168  }
2169  for (Int_t ifidx = fyHiScin[0]; ifidx < (Int_t) fNPaddle[1]; ifidx ++ ){
2170  if ( fScinHitPaddle[1][ifidx] == 1 ){ fHitSweet1Y = -1; }
2171  }
2172 
2173  // *second y plane. first see if there are hits inside the scin region
2174  for (Int_t ifidx = fyLoScin[1]-1; ifidx < fyHiScin[1]; ifidx ++ ){
2175  if ( fScinHitPaddle[3][ifidx] == 1 ){
2176  fHitSweet2Y = 1;
2177  fSweet2YScin = ifidx + 1;
2178  }
2179  }
2180 
2181  // * next make sure nothing fired outside the good region
2182  for (Int_t ifidx = 0; ifidx < fyLoScin[1]-1; ifidx ++ ){
2183  if ( fScinHitPaddle[3][ifidx] == 1 ){ fHitSweet2Y = -1; }
2184  }
2185  for (Int_t ifidx = fyHiScin[1]; ifidx < (Int_t) fNPaddle[3]; ifidx ++ ){
2186  if ( fScinHitPaddle[3][ifidx] == 1 ){ fHitSweet2Y = -1; }
2187  }
2188 
2190 
2191  // * now define a 3/4 or 4/4 trigger of only good scintillators the value
2192  // * is specified in htracking
2194  fGoodScinHits = 1;
2195  for (Int_t ifidx = fxLoScin[0]; ifidx < fxHiScin[0]; ifidx ++ ){
2196  if ( fSweet1XScin == ifidx )
2197  fGoodScinHitsX[ifidx] = 1;
2198  }
2199  }
2200 
2201  // * require front/back hodoscopes be close to each other
2202  if ( ( fGoodScinHits == 1 ) && ( fTrackEffTestNScinPlanes == 4 ) ){
2203  if ( TMath::Abs( fSweet1XScin - fSweet2XScin ) > 3 )
2204  fGoodScinHits = 0;
2205  if ( TMath::Abs( fSweet1YScin - fSweet2YScin ) > 2 )
2206  fGoodScinHits = 0;
2207  }
2208 //
2209 }
2210 //_____________________________________________________________________________
2212 {
2213  Int_t Ntracks = tracks.GetLast()+1; // Number of reconstructed tracks
2214  Double_t hitPos;
2215  Double_t hitDistance;
2216  Int_t ih=0;
2217  for (Int_t itrk=0; itrk<Ntracks; itrk++) {
2218  THaTrack* theTrack = static_cast<THaTrack*>( tracks[itrk] );
2219  if (theTrack->GetIndex()==0) {
2220  fBeta=theTrack->GetBeta();
2221  fFPTimeAll=theTrack->GetFPTime();
2222  Double_t shower_track_enorm = theTrack->GetEnergy()/theTrack->GetP();
2223  for (Int_t ip = 0; ip < fNumPlanesBetaCalc; ip++ ){
2224  Double_t pl_xypos=0;
2225  Double_t pl_zpos=0;
2226  Int_t num_good_pad=0;
2227  TClonesArray* hodoHits = fPlanes[ip]->GetHits();
2228  for (Int_t iphit = 0; iphit < fPlanes[ip]->GetNScinHits(); iphit++ ){
2229  THcHodoHit *hit = fTOFPInfo[ih].hit;
2230  if ( fTOFCalc[ih].good_tdc_pos && fTOFCalc[ih].good_tdc_neg ) {
2231  Bool_t sh_pid=(shower_track_enorm > fTOFCalib_shtrk_lo && shower_track_enorm < fTOFCalib_shtrk_hi);
2233  // cer_pid is true if there is no Cherenkov detector
2235  if(fDumpTOF && Ntracks==1 && fGoodEventTOFCalib && sh_pid && beta_pid && cer_pid) {
2236  fDumpOut << fixed << setprecision(2);
2237  fDumpOut << showpoint << " 1" << setw(3) << ip+1 << setw(3) << hit->GetPaddleNumber() << setw(10) << hit->GetPosTDC()*fScinTdcToTime << setw(10) << fTOFPInfo[ih].pathp << setw(10) << fTOFPInfo[ih].zcor << setw(10) << fTOFPInfo[ih].time_pos << setw(10) << hit->GetPosADC() << endl;
2238  fDumpOut << showpoint << " 2" << setw(3) << ip+1 << setw(3) << hit->GetPaddleNumber() << setw(10) << hit->GetNegTDC()*fScinTdcToTime << setw(10) << fTOFPInfo[ih].pathn << setw(10) << fTOFPInfo[ih].zcor << setw(10) << fTOFPInfo[ih].time_neg << setw(10) << hit->GetNegADC() << endl;
2239  }
2240  Int_t padind = ((THcHodoHit*)hodoHits->At(iphit))->GetPaddleNumber()-1;
2241  pl_xypos+=fPlanes[ip]->GetPosCenter(padind)+ fPlanes[ip]->GetPosOffset();
2242  pl_zpos+=fPlanes[ip]->GetZpos()+ (padind%2)*fPlanes[ip]->GetDzpos();
2243  num_good_pad++;
2244  }
2245  ih++;
2246  // cout << ip << " " << iphit << " " << fGoodFlags[itrk][ip][iphit].goodScinTime << " " << fGoodFlags[itrk][ip][iphit].goodTdcPos << " " << fGoodFlags[itrk][ip][iphit].goodTdcNeg << endl;
2247  }
2248  hitDistance=kBig;
2249  if (num_good_pad !=0 ) {
2250  pl_xypos=pl_xypos/num_good_pad;
2251  pl_zpos=pl_zpos/num_good_pad;
2252  hitPos = theTrack->GetY() + theTrack->GetPhi()*pl_zpos ;
2253  if (ip%2 == 0) hitPos = theTrack->GetX() + theTrack->GetTheta()*pl_zpos ;
2254  hitDistance = hitPos - pl_xypos;
2255  fPlanes[ip]->SetTrackXPosition(theTrack->GetX() + theTrack->GetTheta()*pl_zpos );
2256  fPlanes[ip]->SetTrackYPosition(theTrack->GetY() + theTrack->GetPhi()*pl_zpos );
2257  }
2258  // cout << " ip " << ip << " " << hitPos << " " << pl_xypos << " " << pl_zpos << " " << hitDistance << endl;
2259  fPlanes[ip]->SetHitDistance(hitDistance);
2260  }
2261  if(ih>0&&fDumpTOF && Ntracks==1 && fGoodEventTOFCalib && shower_track_enorm > fTOFCalib_shtrk_lo && shower_track_enorm < fTOFCalib_shtrk_hi ) {
2262  fDumpOut << "0 " << endl;
2263  }
2264  }
2265  } //over tracks
2266  //
2267  return 0;
2268 }
2269 //_____________________________________________________________________________
2271  // GN: Return the index of a scintillator given the plane # and the paddle #
2272  // This assumes that both planes and
2273  // paddles start counting from 0!
2274  // Result also counts from 0.
2275  return fNPlanes*nPaddle+nPlane;
2276 }
2277 //_____________________________________________________________________________
2278 Int_t THcHodoscope::GetScinIndex( Int_t nSide, Int_t nPlane, Int_t nPaddle ) {
2279  return nSide*fMaxHodoScin+fNPlanes*nPaddle+nPlane-1;
2280 }
2281 //_____________________________________________________________________________
2283  return fPathLengthCentral;
2284 }
2285 
2286 //_____________________________________________________________________________
2287 Int_t THcHodoscope::End(THaRunBase* run)
2288 {
2289  MissReport(Form("%s.%s", GetApparatus()->GetName(), GetName()));
2290  return 0;
2291 }
Double_t * fHodoSigmaPos
Definition: THcHodoscope.h:228
Double_t GetPosCorrectedTime() const
Definition: THcHodoHit.h:42
TClonesArray * GetHits()
std::string GetName(const std::string &scope_name)
A single plane of scintillators.
virtual Int_t Fill(Double_t x)
Int_t fHitSweet1Y
Definition: THcHodoscope.h:283
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
Double_t * fHodoNegAdcTimeWindowMin
Definition: THcHodoscope.h:194
Double_t * fHodoPosInvAdcOffset
Definition: THcHodoscope.h:212
void EstimateFocalPlaneTime(void)
void OriginalTrackEffTest(void)
Double_t fTofTolerance
Definition: THcHodoscope.h:185
std::vector< Int_t > fThreeScin
Definition: THcHodoscope.h:392
Double_t GetPosADCtime() const
Definition: THcHodoHit.h:35
Double_t fTOFCalib_cer_lo
Definition: THcHodoscope.h:273
void Setup(const char *name, const char *description)
Bool_t fGoodEventTOFCalib
Definition: THcHodoscope.h:279
Double_t * fPlaneSpacing
Definition: THcHodoscope.h:254
Double_t * fAdcTdcOffset
Definition: THcHodoscope.h:260
Double_t fNCerNPE
Definition: THcHodoscope.h:257
void DeleteArrays()
Int_t fTofUsingInvAdc
Definition: THcHodoscope.h:211
Int_t GetScinIndex(Int_t nPlane, Int_t nPaddle)
std::vector< std::vector< Int_t > > fClustSize
Definition: THcHodoscope.h:386
Int_t GetLast() const
Int_t fTrackBetaIncludeSinglePmtHits
Definition: THcHodoscope.h:192
Int_t fTrackEffTestNScinPlanes
Definition: THcHodoscope.h:263
const char Option_t
Class for gas Cherenkov detectors.
Definition: THcCherenkov.h:16
Double_t fTOFCalib_beta_hi
Definition: THcHodoscope.h:275
Class representing a single hit for the Hodoscopes.
Definition: THcHodoHit.h:16
virtual Int_t End(THaRunBase *run=0)
Double_t fOffsetTime
Definition: THcHodoscope.h:162
void SetClusterUsedFlag(Int_t ic, Double_t flag)
Double_t GetPosTOFCorrectedTime() const
Definition: THcHodoHit.h:44
void MissReport(const char *name)
Definition: THcHitList.cxx:557
Double_t fScinTdcMin
Definition: THcHodoscope.h:189
void SetScinYPos(Double_t f)
Double_t * fHodo_LCoeff
Definition: THcHodoscope.h:222
Double_t * fHodoNegMinPh
Definition: THcHodoscope.h:204
std::vector< std::vector< std::vector< GoodFlags > > > fGoodFlags
Definition: THcHodoscope.h:404
Int_t * fNPlaneTime
Definition: THcHodoscope.h:300
Double_t GetPathLengthCentral()
std::vector< Int_t > fGoodScinHitsX
Definition: THcHodoscope.h:393
Int_t fADC_RefTimeCut
Definition: THcHodoscope.h:146
Int_t GetNegTDC() const
Definition: THcHodoHit.h:41
void SetFpTime(Double_t f)
Double_t * fHodoPosPhcCoeff
Definition: THcHodoscope.h:205
UInt_t GetEvNum() const
Double_t GetPosADC() const
Definition: THcHodoHit.h:31
int Int_t
bool Bool_t
std::vector< std::vector< Double_t > > fdEdX
Definition: THcHodoscope.h:382
std::vector< TOFPInfo > fTOFPInfo
Definition: THcHodoscope.h:361
Double_t fScinTdcMax
Definition: THcHodoscope.h:189
Double_t fStartTimeCenter
Definition: THcHodoscope.h:184
virtual Int_t GetNbinsX() const
STL namespace.
Double_t * fHodoPosMinPh
Definition: THcHodoscope.h:203
virtual Double_t TimeWalkCorrection(const Int_t &paddle, const ESide side)
Int_t fSweet2XScin
Definition: THcHodoscope.h:289
TClonesArray * fRawHitList
Definition: THcHitList.h:51
std::vector< Int_t > fNScinHit
Definition: THcHodoscope.h:383
Double_t * fHodoPosSigma
Definition: THcHodoscope.h:200
THcScintillatorPlane ** fPlanes
Definition: THcHodoscope.h:234
Short_t Abs(Short_t d)
Int_t fSweet2YScin
Definition: THcHodoscope.h:290
Double_t fBetaNominal
Definition: THcHodoscope.h:232
Double_t * fHodoVelLight
Definition: THcHodoscope.h:199
Int_t fCosmicFlag
Definition: THcHodoscope.h:186
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
double beta(double x, double y)
virtual void Reset(Option_t *option="")
Double_t fTimeHist_StartTime_Hits
Definition: THcHodoscope.h:169
std::vector< Int_t > fNCluster
Definition: THcHodoscope.h:388
virtual Int_t DefineVariables(EMode mode=kDefine)
Double_t fStartTimeSlop
Definition: THcHodoscope.h:184
Double_t fBetaNoTrk
Definition: THcHodoscope.h:177
void SetTwoGoodTimes(Bool_t flag)
Definition: THcHodoHit.h:64
Double_t fTOFCalib_shtrk_lo
Definition: THcHodoscope.h:271
void TrackEffTest(void)
virtual EStatus Init(const TDatime &run_time)
UInt_t fMaxHodoScin
Definition: THcHodoscope.h:183
std::vector< std::vector< Int_t > > fScinHitPaddle
Definition: THcHodoscope.h:384
Double_t fTdc_Thrs
Definition: THcHodoscope.h:227
Int_t fAnalyzePedestals
Definition: THcHodoscope.h:148
Double_t GetNegADCCorrtime() const
Definition: THcHodoHit.h:38
Double_t GetPosADCpeak() const
Definition: THcHodoHit.h:33
Bool_t GetTwoGoodTimes() const
Definition: THcHodoHit.h:47
Int_t GetPaddleNumber() const
Definition: THcHodoHit.h:49
Double_t fPathLengthCentral
Definition: THcHodoscope.h:188
char ** fPlaneNames
Definition: THcHodoscope.h:190
std::vector< TOFCalc > fTOFCalc
Definition: THcHodoscope.h:378
double pow(double, double)
Int_t fGoodScinHits
Definition: THcHodoscope.h:264
tuple app
virtual Int_t CoarseProcess(TClonesArray &tracks)
virtual void Clear(Option_t *opt="")
const Double_t sigma
Double_t DetermineTimePeak(Int_t FillFlag)
Double_t * fHodoNegAdcTimeWindowMax
Definition: THcHodoscope.h:195
Double_t fTimeHist_StartTime_Peak
Definition: THcHodoscope.h:167
std::vector< std::vector< Double_t > > fClusterXPos
Definition: THcHodoscope.h:390
Int_t * fyHiScin
Definition: THcHodoscope.h:268
void Error(const char *location, const char *msgfmt,...)
virtual Int_t AccumulatePedestals(TClonesArray *rawhits, Int_t nexthit)
Double_t * fHodoPosAdcTimeWindowMax
Definition: THcHodoscope.h:197
Double_t fStartTime
Definition: THcHodoscope.h:160
Int_t GetPosTDC() const
Definition: THcHodoHit.h:40
Double_t GetCerNPE()
void InitHitList(THaDetMap *detmap, const char *hitclass, Int_t maxhits, Int_t tdcref_cut=0, Int_t adcref_cut=0)
Save the electronics module to detector mapping and initialize a hit array of hits of class hitclass...
Definition: THcHitList.cxx:65
Double_t fTimeHist_FpTime_Sigma
Definition: THcHodoscope.h:172
Int_t fNumPlanesBetaCalc
Definition: THcHodoscope.h:187
void SetClusterFlag(Int_t ic, Double_t flag)
Int_t fHitSweet2Y
Definition: THcHodoscope.h:285
Double_t * fHodoNeg_c1
Definition: THcHodoscope.h:224
virtual void SetRange(Int_t first=0, Int_t last=0)
Double_t trackeff_scint_xdiff_max
Definition: THcHodoscope.h:156
virtual Double_t GetMean(Int_t axis=1) const
Double_t GetNegCorrectedTime() const
Definition: THcHodoHit.h:43
virtual Double_t Integral(Option_t *option="") const
Int_t fCheckEvent
Definition: THcHodoscope.h:242
Double_t fTimeHist_FpTime_Peak
Definition: THcHodoscope.h:171
void SetClusterSize(Int_t ic, Double_t size)
UInt_t fMaxScinPerPlane
Definition: THcHodoscope.h:183
Bool_t * fPresentP
Definition: THcHodoscope.h:165
Int_t fSweet1XScin
Definition: THcHodoscope.h:287
Int_t fEventType
Definition: THcHodoscope.h:243
Int_t * fNScinHits
Definition: THcHodoscope.h:299
Double_t GetParticleMass() const
Double_t * fHodoSigmaNeg
Definition: THcHodoscope.h:229
void CalcCluster(void)
unsigned int UInt_t
Int_t * fyLoScin
Definition: THcHodoscope.h:267
std::vector< std::vector< Double_t > > fClustPos
Definition: THcHodoscope.h:387
char * Form(const char *fmt,...)
void CreateMissReportParms(const char *prefix)
Definition: THcHitList.cxx:544
ofstream fDumpOut
Definition: THcHodoscope.h:277
Double_t GetNegADCpeak() const
Definition: THcHodoHit.h:34
Double_t * fHodoPos_c1
Definition: THcHodoscope.h:223
Double_t * fHodoPosInvAdcLinear
Definition: THcHodoscope.h:214
Int_t fdebugprintscinraw
Definition: THcHodoscope.h:261
Double_t * fHodoPos_c2
Definition: THcHodoscope.h:225
Double_t GetRMS(Int_t axis=1) const
Double_t trackeff_scint_ydiff_max
Definition: THcHodoscope.h:155
virtual void Clear(Option_t *opt="")
std::vector< std::vector< Int_t > > fClusterSize
Definition: THcHodoscope.h:389
tuple list
Definition: SConscript.py:9
void SetNumberClusters(Int_t nclus)
Int_t fSweet1YScin
Definition: THcHodoscope.h:288
const Bool_t kFALSE
virtual Int_t ReadDatabase(const TDatime &date)
Double_t fScinTdcToTime
Definition: THcHodoscope.h:184
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition: THcGlobals.h:12
UInt_t fTotHodScin
Definition: THcHodoscope.h:183
Double_t * fHodoSlop
Definition: THcHodoscope.h:258
Double_t fTOFCalib_beta_lo
Definition: THcHodoscope.h:274
Double_t fTimeHist_StartTime_Sigma
Definition: THcHodoscope.h:168
virtual Int_t Decode(const THaEvData &)
Double_t fTimeHist_StartTime_NumPeaks
Definition: THcHodoscope.h:166
Int_t * fTdcOffset
Definition: THcHodoscope.h:259
Double_t * fHodoPosAdcTimeWindowMin
Definition: THcHodoscope.h:196
Double_t * fSumPlaneTime
Definition: THcHodoscope.h:297
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
double Double_t
Int_t fNHodoscopes
Definition: THcHodoscope.h:269
Double_t GetPosADCCorrtime() const
Definition: THcHodoHit.h:37
Int_t fHitSweet1X
Definition: THcHodoscope.h:282
virtual Int_t ProcessHits(TClonesArray *rawhits, Int_t nexthit)
Double_t fPartMass
Definition: THcHodoscope.h:231
Double_t * fFPTime
Definition: THcHodoscope.h:294
Int_t fDebug
Definition: THcParmList.cxx:62
Double_t * fHodoPosInvAdcAdc
Definition: THcHodoscope.h:216
void SetTrackXPosition(Double_t f)
Double_t fADCStartTime
Definition: THcHodoscope.h:161
Double_t * fHodoNegInvAdcLinear
Definition: THcHodoscope.h:215
Double_t * fHodoNegSigma
Definition: THcHodoscope.h:201
UInt_t GetEvType() const
Double_t GetBetaAtPcentral() const
Double_t * fHodoNegPhcCoeff
Definition: THcHodoscope.h:206
Bool_t GetHasCorrectedTimes() const
Definition: THcHodoHit.h:48
ClassImp(THcDCLookupTTDConv) THcDCLookupTTDConv
void SetCluster(Int_t ic, Double_t pos)
UInt_t fNRawHits
Definition: THcHitList.h:45
Bool_t fGoodStartTime
Definition: THcHodoscope.h:159
virtual Int_t ApplyCorrections(void)
Double_t GetNegADCtime() const
Definition: THcHodoHit.h:36
Double_t * fHodoNeg_c2
Definition: THcHodoscope.h:226
Double_t GetNegADC() const
Definition: THcHodoHit.h:32
UInt_t * fNPaddle
Definition: THcHodoscope.h:191
void SetScinXPos(Double_t f)
Double_t * fHodoCableFit
Definition: THcHodoscope.h:221
std::vector< Int_t > fNClust
Definition: THcHodoscope.h:385
virtual ~THcHodoscope()
Int_t fHitSweet2X
Definition: THcHodoscope.h:284
Double_t fBetaNoTrkChiSq
Definition: THcHodoscope.h:178
Double_t * fHodoPosTimeOffset
Definition: THcHodoscope.h:207
Short_t Max(Short_t a, Short_t b)
Double_t * fHodoNegInvAdcAdc
Definition: THcHodoscope.h:217
Double_t GetScinCorrectedTime() const
Definition: THcHodoHit.h:46
void SetHitDistance(Double_t f)
Int_t * fxHiScin
Definition: THcHodoscope.h:266
void SetTrackYPosition(Double_t f)
virtual Int_t FineProcess(TClonesArray &tracks)
Double_t fTimeHist_FpTime_NumPeaks
Definition: THcHodoscope.h:170
Double_t fTimeHist_FpTime_Hits
Definition: THcHodoscope.h:173
Double_t fTOFCalib_shtrk_hi
Definition: THcHodoscope.h:272
Double_t fRatio_xpfp_to_xfp
Definition: THcHodoscope.h:154
Int_t fTDC_RefTimeCut
Definition: THcHodoscope.h:145
R__EXTERN class THcParmList * gHcParms
Definition: THcGlobals.h:11
Double_t * fPlaneCenter
Definition: THcHodoscope.h:253
Double_t Sqrt(Double_t x)
TObject * At(Int_t idx) const
static const Double_t kBig
Definition: THcFormula.cxx:31
Int_t * fHodoPosPedLimit
Definition: THcHodoscope.h:209
Int_t * fHodoNegPedLimit
Definition: THcHodoscope.h:210
Double_t * fHodoNegTimeOffset
Definition: THcHodoscope.h:208
const Bool_t kTRUE
Generic hodoscope consisting of multiple planes with multiple paddles with phototubes on both ends...
Definition: THcHodoscope.h:37
Double_t * fHodoVelFit
Definition: THcHodoscope.h:220
string fTOFDumpFile
Definition: THcHodoscope.h:278
Double_t GetPosCenter(Int_t PaddleNo)
Double_t fNormETot
Definition: THcHodoscope.h:256
Double_t GetNegTOFCorrectedTime() const
Definition: THcHodoHit.h:45
Double_t fBeta
Definition: THcHodoscope.h:175
Double_t fFPTimeAll
Definition: THcHodoscope.h:163
Double_t * fHodoNegInvAdcOffset
Definition: THcHodoscope.h:213
virtual Int_t DecodeToHitList(const THaEvData &evdata, Bool_t suppress=kFALSE)
Populate the hitlist from the raw event data.
Definition: THcHitList.cxx:191
TAxis * GetXaxis()
THcCherenkov * fCherenkov
Definition: THcHodoscope.h:143
std::vector< std::vector< Double_t > > fClusterYPos
Definition: THcHodoscope.h:391
Bool_t * fGoodPlaneTime
Definition: THcHodoscope.h:302
Int_t * fxLoScin
Definition: THcHodoscope.h:265
A standard Hall C spectrometer apparatus.