/search.css" rel="stylesheet" type="text/css"/> /search.js">
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

Public Member Functions | Static Public Member Functions | Private Attributes
TimeRecoTool Class Reference

#include <TimeRecoTool.h>

Inheritance diagram for TimeRecoTool:
Inheritance graph
[legend]
Collaboration diagram for TimeRecoTool:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 TimeRecoTool (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~TimeRecoTool ()
virtual StatusCode reconstruct (const DayaBay::CalibReadout &readout, DayaBay::RecTrigger &recTrigger)
virtual StatusCode initialize ()
virtual StatusCode finalize ()

Static Public Member Functions

static const InterfaceID & interfaceID ()

Private Attributes

std::string m_templateFileName
std::string m_cableSvcName_t
ICableSvcm_cableSvc_t
std::string m_calibDataSvcName_t
IPmtCalibSvcm_calibDataSvc_t
std::string m_pmtGeomSvcName_t
IPmtGeomInfoSvcm_pmtGeomSvc_t
std::string m_ChannelQualitySvcName_t
IChannelQualitySvcm_chanqualSvc_t
bool m_useChannelQualitySvc_t
std::string m_AdSimpleLocation_t
int nocqCount
TFile * savefile
TCanvas * ccharge
TCanvas * ctime
TH1F * hchantime1
TH1F * hchantime0
TH1F * hchancharge1
TH1F * hchancharge0
TF1 * fittime1
TF1 * fittime0
TF1 * reftime0
TF1 * reftime1
TH1F * hsteps
TH2F * hstepsVdDis
TH2F * hstepsVchi2
TH1F * hdifftime0
bool Debug_t
bool m_debugOutput_t

Detailed Description

Definition at line 33 of file TimeRecoTool.h.


Constructor & Destructor Documentation

TimeRecoTool::TimeRecoTool ( const std::string &  type,
const std::string &  name,
const IInterface *  parent 
)

Definition at line 84 of file TimeRecoTool.cc.

  : GaudiTool(type,name,parent)
  , m_cableSvc_t(0)
  , m_calibDataSvc_t(0)
  , m_pmtGeomSvc_t(0)
  , m_chanqualSvc_t(0)
{


    declareInterface< IReconTool >(this) ;   
    
    declareProperty("CableSvcName",m_cableSvcName_t="CableSvc",
                    "Name of service to map between detector, hardware, and electronic IDs");
    declareProperty("CalibDataSvcName", m_calibDataSvcName_t="DybPmtCalibSvc",
                    "Name of calibration data service");
    declareProperty("PmtGeomSvcName", m_pmtGeomSvcName_t="PmtGeomInfoSvc",
                    "Name of PMT geometry data service");
    declareProperty("DebugOutput", m_debugOutput_t=false,
                    "Output debug informaiton");
    declareProperty("ChannelQualitySvcName", m_ChannelQualitySvcName_t="DybChannelQualitySvc",
                    "Name of channel quality service");
    declareProperty("UseChannelQualitySvc", m_useChannelQualitySvc_t=true,
                    "Use ChannelQualityService to check if the PMTs are alive. Turn it off if it's too slow");

}
TimeRecoTool::~TimeRecoTool ( ) [virtual]

Definition at line 112 of file TimeRecoTool.cc.

{}

Member Function Documentation

StatusCode TimeRecoTool::reconstruct ( const DayaBay::CalibReadout readout,
DayaBay::RecTrigger recTrigger 
) [virtual]

Implements IReconTool.

Definition at line 197 of file TimeRecoTool.cc.

{
  
  entry ++;
  
  if( !readout.detector().isAD() ){
    debug() << "Not an AD readout; ignoring detector " 
            << readout.detector().detName() << endreq;
    recTrigger.setPositionStatus( ReconStatus::kNotProcessed );
    return StatusCode::SUCCESS;
  }

  const CalibReadoutPmtCrate* pmtReadout 
    = dynamic_cast<const CalibReadoutPmtCrate*>(&readout);
  if(!pmtReadout){
    error() << "Incorrect type of readout crate for detector "
            << readout.detector().detName() << endreq;
    recTrigger.setPositionStatus( ReconStatus::kBadReadout );
    return StatusCode::FAILURE;
  }


  // Context for this data
  int task = 0;
  ServiceMode svcMode(readout.header()->context(), task);

  Int_t nActivePMTs_t = 0;

  if (m_useChannelQualitySvc_t){
    const IChannelQuality * cq = m_chanqualSvc_t->channelQuality(svcMode);
    if (!cq){
      nocqCount ++;
      if (nocqCount <= NOCQ_ERRORS){
        warning() << "\tGot null channel quality object for " << readout.header()->context().AsString()
                  << ". Assuming all channels are good." << endreq;
      }
      if (nocqCount == NOCQ_ERRORS){
        warning() << "\tFurther warning for channel quality service will be suppressed" << endreq;
      }

      nActivePMTs_t = 192;
      for (Int_t iring = 0; iring < 8; iring++){
        for (Int_t icolumn = 0; icolumn < 24; icolumn++){
          pmt_flag_t[iring][icolumn] = 1;
        }
      }

    }else{
      
      // loop over channels and check if PMTs are active
      IChannelQuality::ChannelSet_t chans = cq->channels();
      for (int i = 0; i < chans.size(); i++){
        DayaBay::FeeChannelId chid = chans[i];
        AdPmtSensor pmtId = m_cableSvc_t->adPmtSensor(chid, svcMode);
        if (pmtId.is8inch()){
          if (cq->good(chid)){
            pmt_flag_t[pmtId.ring()-1][pmtId.column()-1] = 1;
            nActivePMTs_t ++;
          }else{
            pmt_flag_t[pmtId.ring()-1][pmtId.column()-1] = -1;
            //    std::cout << "\tGod bad one: ring " << pmtId.ring() << " column " << pmtId.column() << std::endl;
          
          }
        }
      }
    }
    //    std::cout << "Number of active PMTs = " << nActivePMTs << std::endl;
  }else{ // all PMTs are assumed to be active
    nActivePMTs_t = 192;
    for (Int_t iring = 0; iring < 8; iring++){
      for (Int_t icolumn = 0; icolumn < 24; icolumn++){
        pmt_flag_t[iring][icolumn] = 1;
      }
    }
  }
 
  for (Int_t iring = 0; iring < 8; iring++){
    for (Int_t icolumn = 0; icolumn < 24; icolumn++){ 
      
      pmt_charge_t[iring][icolumn] = 0;
      pmt_time_t[iring][icolumn] = 0;
      pmt_error_t[iring][icolumn] = 0; 
      deltime_t[iring][icolumn] = 0;
    
    }
  }
  
  //Loop over channels 
  //******************
  Int_t nHits_t = 0;
  total_pe_t = 0;
  CalibReadoutPmtCrate::PmtChannelReadouts::const_iterator chanIter, 
    chanEnd = pmtReadout->channelReadout().end();
  for(chanIter=pmtReadout->channelReadout().begin(); chanIter != chanEnd; 
      chanIter++){
    const CalibReadoutPmtChannel& channel = *chanIter;
    AdPmtSensor pmtId(channel.pmtSensorId().fullPackedData());
    if( !pmtId.is8inch() ){ continue; }  // Ignore non 8-inch pmt channels
    const PmtCalibData* pmtCalib = m_calibDataSvc_t->pmtCalibData(pmtId, svcMode);
    if( !pmtCalib ){ 
      error() << "No calibration data for pmt ID: " << pmtId << endreq;
      recTrigger.setPositionStatus( ReconStatus::kBadReadout );
      return StatusCode::FAILURE;
    }
    const CLHEP::Hep3Vector& pmtPos 
      = m_pmtGeomSvc_t->get(pmtId.fullPackedData())->localPosition();
    double peakAdc=channel.earliestCharge(-1650,-1250);
    double firstTdc = channel.earliestTime();
    if(pmt_flag_t[pmtId.ring()-1][pmtId.column()-1] == 1){

      pmt_charge_t[pmtId.ring()-1][pmtId.column()-1] = peakAdc;
      pmt_time_t[pmtId.ring()-1][pmtId.column()-1] = firstTdc;

      pmt_x_t[pmtId.ring()-1][pmtId.column()-1] = pmtPos[0];
      pmt_y_t[pmtId.ring()-1][pmtId.column()-1] = pmtPos[1];
      pmt_z_t[pmtId.ring()-1][pmtId.column()-1] = pmtPos[2];

      if (peakAdc>0){
        nHits_t++;
        total_pe_t += peakAdc;
      }
 
    }
  }
  
  if(nHits_t<1){
    CLHEP::HepLorentzVector position(1e10, 1e10, 1e10, 1e10); // set dummy number in case someone try to use it
    recTrigger.setPosition( position );
    recTrigger.setPositionQuality( 1e10 );
    recTrigger.setPositionStatus( ReconStatus::kNoHits ); //No Hits
    return StatusCode::SUCCESS;
  }

 // correction of the total pe based on the number of active PMTs
  total_pe_t *= 192./ (Double_t)nActivePMTs_t;

  //Initialize final position
  final_pos_t[0] = 1e10;
  final_pos_t[1] = 1e10;
  final_pos_t[2] = 1e10;
  final_pos_t[3] = 1e10;

  //Load Prefit Vertex (AdSimple), otherwise use a simple Charge guess
  //****************************************************************** 
  if(recTrigger.positionStatus()==ReconStatus::kGood){
    const CLHEP::HepLorentzVector oldpos = recTrigger.position();
    initial_pos_t[0] = oldpos.x();
    initial_pos_t[1] = oldpos.y();
    initial_pos_t[2] = oldpos.z();
    initial_pos_t[3] = oldpos.t();

    //Check if initial r is too close to boundary
    Double_t r = sqrt(pow(initial_pos_t[0],2) + pow(initial_pos_t[1],2));
    if (r > 1800){
      Double_t phi = atan2(initial_pos_t[1],initial_pos_t[0]);
      initial_pos_t[0] = 1800*cos(phi);
      initial_pos_t[1] = 1800*sin(phi);
    }
    
    //Check if initial z is too close to boundary
    if (initial_pos_t[2] > 1700) initial_pos_t[2] =  1700;
    else if (initial_pos_t[2] < -1700) initial_pos_t[2] = -1700;
    
    initial_energy_t = recTrigger.energy();
    if (initial_energy_t<0) initial_energy_t = total_pe_t / 163;

  }
  else{
    SimpleStart();
    initial_energy_t = total_pe_t / 163;
  }

  Double_t minPE;
  if (initial_energy_t > 2.5) minPE = 0.5;
  else minPE = initial_energy_t/5;
  
  Int_t nGoodHits_t = 0;
  for (Int_t iring = 0; iring < 8; iring++){
    for (Int_t icolumn = 0; icolumn < 24; icolumn++){
      if (pmt_flag_t[iring][icolumn] == 1){

        //Eliminate low charge and nohit PMT
        if (pmt_charge_t[iring][icolumn] < minPE) 
          pmt_flag_t[iring][icolumn] = -1;
        
        //Set PMT timing error for good PMT     
        else{
          nGoodHits_t++;
          if (pmt_charge_t[iring][icolumn]<41.08)
            pmt_error_t[iring][icolumn] =  (0.43 + 31.5 * pow(pmt_charge_t[iring][icolumn]*18.5,-0.67)) / sqrt(5);
          else
            pmt_error_t[iring][icolumn] =  0.80 / sqrt(5);  
        }
      
      }
    }
  }
 
  if(nGoodHits_t<3){
    recTrigger.setPositionStatus( ReconStatus::kNotConverged ); //Needs at least 3 good PMT
    return StatusCode::SUCCESS;
  }  
  
  //Calculate deltime and select PMTs
  //*********************************
  CalculateDeltime(initial_pos_t[0],initial_pos_t[1],initial_pos_t[2]);
  SelectPMT();
  
  Int_t iSelect = 0;
  Bool_t SuccessFlag = 0;
  Int_t nFitPMT_t = 0;

  while (SuccessFlag == 0 && iSelect<maxselect_t){
    
    if(nPMT_t[iSelect] > 0){
      nFitPMT_t += nPMT_t[iSelect];
      
      if(nFitPMT_t > 2) //Do fit only if there are at least 3 PMTs
        IterativeFit(nFitPMT_t,maxstep_t,initial_pos_t[0],initial_pos_t[1],initial_pos_t[2], iSelect);
    }

    if (sqrt(pow(final_pos_t[0],2)+pow(final_pos_t[1],2)) < 2300 && final_pos_t[2] < 2200 && final_pos_t[2] > -2200) SuccessFlag = 1; //Good convergence
    else iSelect++;
    
  }
  
  if(SuccessFlag == 0){ //Check if convergence occurred outsite physical boundary
    if (sqrt(pow(final_pos_t[0],2)+pow(final_pos_t[1],2)) < 10*r_phys_max_t && final_pos_t[2] < 10*z_phys_max_t && final_pos_t[2] > 10*z_phys_min_t)
      SuccessFlag = 1; //Convergenced
  }

  if(SuccessFlag == 0){  //Did not converge
    recTrigger.setPositionStatus( ReconStatus::kNotConverged );
    return StatusCode::SUCCESS;
  }
 
  CLHEP::HepLorentzVector position(final_pos_t[0], final_pos_t[1], final_pos_t[2], final_pos_t[3]);
  recTrigger.setPosition( position );
  recTrigger.setPositionQuality( chi2_sum_t(final_pos_t[0],final_pos_t[1],final_pos_t[2],final_pos_t[3],iSelect) );
  recTrigger.setPositionStatus( ReconStatus::kGood );
  
  //Debugging
  //*********

  if (Debug_t){
    if ((readout.triggerType() & DayaBay::Trigger::kCalib)==268435472){
      hsteps->Fill(totalsteps);
      hstepsVdDis->Fill(totalsteps,sqrt((final_pos_t[0]-initial_pos_t[0])*(final_pos_t[0]-initial_pos_t[0])+(final_pos_t[1]-initial_pos_t[1])*(final_pos_t[1]-initial_pos_t[1])+(final_pos_t[2]-initial_pos_t[2])*(final_pos_t[2]-initial_pos_t[2])));
      hstepsVchi2->Fill(totalsteps,recTrigger.positionQuality());
    }

    if ((readout.triggerType() & DayaBay::Trigger::kCalib)==268435472){ // && recTrigger.positionQuality()>0.3){ // && recTrigger.positionStatus() > 2){
      hchantime1->Reset();
      hchantime0->Reset();
      hchancharge1->Reset();
      hchancharge0->Reset();
      
      hdifftime0->Reset();
      
      for (Int_t iring = 0; iring < 8; iring++){
        for (Int_t icolumn = 0; icolumn < 24; icolumn++){
          Int_t ichan = (icolumn+1)+ 24*iring; 
          
          if(pmt_flag_t[iring][icolumn] == 1){
            hchantime1->SetBinContent(ichan,pmt_time_t[iring][icolumn] - median_deltime_t);  
            hchantime1->SetBinError(ichan,pmt_error_t[iring][icolumn]);
            hchancharge1->SetBinContent(ichan,pmt_charge_t[iring][icolumn]);  
          }
          else if (pmt_flag_t[iring][icolumn] == 0){
            hchantime0->SetBinContent(ichan,pmt_time_t[iring][icolumn] - median_deltime_t);  
            hchantime0->SetBinError(ichan,pmt_error_t[iring][icolumn]);
            hchancharge0->SetBinContent(ichan,pmt_charge_t[iring][icolumn]); 
          }
          
          if(pmt_charge_t[iring][icolumn]>0.2) hdifftime0->Fill(deltime_t[iring][icolumn] - median_deltime_t);
        }
      }
      
      fittime1->FixParameter(0,final_pos_t[3] - median_deltime_t); 
      fittime1->FixParameter(1,final_pos_t[0]);
      fittime1->FixParameter(2,final_pos_t[1]);  
      fittime1->FixParameter(3,final_pos_t[2]);
      
      fittime1->SetLineWidth(1);
      fittime1->SetLineColor(2);
      hchantime1->Fit(fittime1,"Q");
      
      fittime0->FixParameter(0,final_pos_t[3] - median_deltime_t ); 
      fittime0->FixParameter(1,initial_pos_t[0]);
      fittime0->FixParameter(2,initial_pos_t[1]);
      fittime0->FixParameter(3,initial_pos_t[2]);
      
      fittime0->SetLineWidth(1);
      fittime0->SetLineColor(4);
      fittime0->SetLineStyle(9);
      hchantime1->Fit(fittime0,"+Q");
      
      reftime0->FixParameter(0,early_deltime_t - median_deltime_t);
      reftime0->SetLineWidth(1);
      reftime0->SetLineColor(2);
      reftime0->SetLineStyle(2);
      hchantime1->Fit(reftime0,"+Q");
      
      hchantime1->SetTitle(DayaBay::Trigger::AsString(readout.triggerType()));
      
      Char_t hname[1024];
      sprintf(hname, "%s_%i_hchantime1", DayaBay::Trigger::AsString(readout.triggerType()), entry);
      hchantime1->SetName(hname);
      ctime->SetName(hname);
      
      Char_t hname2[1024];
      sprintf(hname2, "%s_%i_hchancharge1", DayaBay::Trigger::AsString(readout.triggerType()), entry);
      hchancharge1->SetName(hname2);
      ccharge->SetName(hname2);
      
      Char_t hname3[1024];
      sprintf(hname3, "%s_%i_hdifftime0", DayaBay::Trigger::AsString(readout.triggerType()), entry);
      hdifftime0->SetName(hname3);
      
      Char_t ytitle[1024];
      sprintf(ytitle, "Initial t=%.2f,x=%.2f,y=%.2f,z=%.2f ", 0.0, initial_pos_t[0] ,initial_pos_t[1] ,  initial_pos_t[2] );
      hchantime1->GetYaxis()->SetTitle(ytitle);
      hchancharge1->GetYaxis()->SetTitle(ytitle);
      
      Char_t xtitle[1024];
      sprintf(xtitle, "(t=%.2f,x=%.2f,y=%.2f,z=%.2f)", final_pos_t[3] - median_deltime_t, final_pos_t[0], final_pos_t[1], final_pos_t[2]);
      hchantime1->GetXaxis()->SetTitle(xtitle);
      hchancharge1->GetXaxis()->SetTitle(xtitle);
      
      Char_t early_xtitle[1024];
      sprintf(early_xtitle, "(early_deltime_t=%.2f)", early_deltime_t - median_deltime_t);
      hdifftime0->GetXaxis()->SetTitle(early_xtitle);
      
      TDirectory * dir = gDirectory;
      
      savefile->cd();
      
      ctime->cd();
      hchantime1->Draw();
      hchantime0->Draw("same");
      ctime->Write();
      
      ccharge->cd();
      hchancharge1->Draw();
      hchancharge0->Draw("same");
      //htemplatecharge1->Draw("same");
      //htemplatecharge0->Draw("same");
      ccharge->Write();
      
      hdifftime0->Write();
      
      dir->cd();
    }
  }
  
  return StatusCode::SUCCESS;
}
StatusCode TimeRecoTool::initialize ( ) [virtual]

Definition at line 114 of file TimeRecoTool.cc.

{
  
  entry = 0;

  // Get Cable Service
  m_cableSvc_t = svc<ICableSvc>(m_cableSvcName_t,true);
  m_calibDataSvc_t = svc<IPmtCalibSvc>(m_calibDataSvcName_t,true);
  m_pmtGeomSvc_t = svc<IPmtGeomInfoSvc>(m_pmtGeomSvcName_t,true);

  // channel quality service
  m_chanqualSvc_t = svc<IChannelQualitySvc>(m_ChannelQualitySvcName_t,true);
  nocqCount=0;

  Debug_t = m_debugOutput_t;
  
  if (Debug_t) {
    
    //Debugging tools
    TDirectory * dir = gDirectory;
    savefile = new TFile("debug.root","RECREATE");

    //fittime1 = new TF1("fittime1","([0]+(1.465/300)*sqrt(pow(([1]*cos([2])+2250*cos(((x-.5)/12)*3.141592)),2)+pow(([1]*sin([2])+2250*sin(((x-.5)/12)*3.141592)),2)+pow(([3]-500.0*(int((x-1)/24)-3.5)),2)))",1,192); 
    //fittime0 = new TF1("fittime0","([0]+(1.465/300)*sqrt(pow(([1]*cos([2])+2250*cos(((x-.5)/12)*3.141592)),2)+pow(([1]*sin([2])+2250*sin(((x-.5)/12)*3.141592)),2)+pow(([3]-500.0*(int((x-1)/24)-3.5)),2)))",1,192); 

    fittime1 = new TF1("fittime1","([0]+(1.465/300)*sqrt(pow(([1]+2250*cos(((x-.5)/12)*3.141592)),2)+pow(([2]+2250*sin(((x-.5)/12)*3.141592)),2)+pow(([3]-500.0*(int((x-1)/24)-3.5)),2)))",1,192); 
    fittime0 = new TF1("fittime0","([0]+(1.465/300)*sqrt(pow(([1]+2250*cos(((x-.5)/12)*3.141592)),2)+pow(([2]+2250*sin(((x-.5)/12)*3.141592)),2)+pow(([3]-500.0*(int((x-1)/24)-3.5)),2)))",1,192); 
    reftime0 = new TF1("reftime0","[0]",1,192); 
    reftime1 = new TF1("reftime1","[0]",1,192); 

    hchantime1 = new TH1F("hchantime1","earlytime vs channum",192,0.5,192.5);
    hchantime0 = new TH1F("hchantime0","earlytime vs channum",192,0.5,192.5);
    hchancharge1 = new TH1F("hchancharge1","charge vs channum",192,0,192);
    hchancharge0 = new TH1F("hchancharge0","charge vs channum",192,0,192);
    
    hdifftime0 = new TH1F("hdifftime0","hdifftime0",300,-50,250);
    hchantime0->SetLineColor(3);
    hchantime0->SetMarkerStyle(4);

    hchancharge0->SetLineColor(3);

    ccharge = new TCanvas();
    ccharge->SetBorderMode(0);
    ccharge->SetFrameBorderMode(0);
    ccharge->SetFillColor(0);

    ctime = new TCanvas();
    ctime->SetBorderMode(0);
    ctime->SetFrameBorderMode(0);
    ctime->SetFillColor(0);
    
    hsteps = new TH1F("hsteps","hsteps",maxstep_t+1,0,maxstep_t+1);
    hstepsVdDis = new TH2F("hstepsVdDis","hstepsVdDis",maxstep_t+1,0,maxstep_t+1,100,0,1000);
    hstepsVchi2 = new TH2F("hstepsVchi2","hstepsVchi2",maxstep_t+1,0,maxstep_t+1,50,0,1);

    dir->cd();
    
  }
 
  return StatusCode::SUCCESS;
}
StatusCode TimeRecoTool::finalize ( ) [virtual]

Definition at line 176 of file TimeRecoTool.cc.

{
  //Status report for null Channel Quality Service
  if (nocqCount > 0)
    warning() << "******* Got null Channel Quality objects for "<< nocqCount << " times  *********" << endreq;
 
 if (Debug_t){
     
   TDirectory * dir = gDirectory;
   savefile->cd();
   hsteps->Write();
   hstepsVdDis->Write();
   hstepsVchi2->Write();
   dir->cd();

    savefile->Close();
  }

  return StatusCode::SUCCESS;
}

Member Data Documentation

std::string TimeRecoTool::m_templateFileName [private]

Definition at line 51 of file TimeRecoTool.h.

std::string TimeRecoTool::m_cableSvcName_t [private]

Definition at line 54 of file TimeRecoTool.h.

Definition at line 56 of file TimeRecoTool.h.

std::string TimeRecoTool::m_calibDataSvcName_t [private]

Definition at line 59 of file TimeRecoTool.h.

Definition at line 61 of file TimeRecoTool.h.

std::string TimeRecoTool::m_pmtGeomSvcName_t [private]

Definition at line 64 of file TimeRecoTool.h.

Definition at line 66 of file TimeRecoTool.h.

Definition at line 69 of file TimeRecoTool.h.

Definition at line 71 of file TimeRecoTool.h.

Definition at line 72 of file TimeRecoTool.h.

std::string TimeRecoTool::m_AdSimpleLocation_t [private]

Definition at line 75 of file TimeRecoTool.h.

int TimeRecoTool::nocqCount [private]

Definition at line 78 of file TimeRecoTool.h.

TFile* TimeRecoTool::savefile [private]

Definition at line 81 of file TimeRecoTool.h.

TCanvas* TimeRecoTool::ccharge [private]

Definition at line 82 of file TimeRecoTool.h.

TCanvas* TimeRecoTool::ctime [private]

Definition at line 83 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hchantime1 [private]

Definition at line 84 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hchantime0 [private]

Definition at line 85 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hchancharge1 [private]

Definition at line 86 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hchancharge0 [private]

Definition at line 87 of file TimeRecoTool.h.

TF1* TimeRecoTool::fittime1 [private]

Definition at line 88 of file TimeRecoTool.h.

TF1* TimeRecoTool::fittime0 [private]

Definition at line 89 of file TimeRecoTool.h.

TF1* TimeRecoTool::reftime0 [private]

Definition at line 90 of file TimeRecoTool.h.

TF1* TimeRecoTool::reftime1 [private]

Definition at line 91 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hsteps [private]

Definition at line 93 of file TimeRecoTool.h.

TH2F* TimeRecoTool::hstepsVdDis [private]

Definition at line 94 of file TimeRecoTool.h.

TH2F* TimeRecoTool::hstepsVchi2 [private]

Definition at line 95 of file TimeRecoTool.h.

TH1F* TimeRecoTool::hdifftime0 [private]

Definition at line 98 of file TimeRecoTool.h.

bool TimeRecoTool::Debug_t [private]

Definition at line 100 of file TimeRecoTool.h.

Definition at line 102 of file TimeRecoTool.h.


The documentation for this class was generated from the following files:
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Fri May 16 2014 10:12:10 for TimeReco by doxygen 1.7.4