BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//----------------------------------------------------------------//
//                                                                //
//	BrTofRdoModule                                            // 
//                                                                // 
//      Module class that takes digitized data (BrTofDig)         //
//      and calibration parameters  to reconstruct TOF hits       //
//                                                                // 
//                                                                // 
//----------------------------------------------------------------//

//
// $Id: BrTofRdoModule.cxx,v 1.20 2002/08/06 19:07:34 hagel Exp $
// $Author: hagel $
// $Date: 2002/08/06 19:07:34 $ 
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> 
//
//----------------------------------------------------------------//
//

#ifndef WIN32
#include <iostream>
#include <iomanip>
#else
#include <iostream.h>
#include <iomanip.h>
#endif

#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TNtuple
#include "TNtuple.h"
#endif
#ifndef BRAT_BrTofRdoModule
#include "BrTofRdoModule.h"
#endif
#ifndef BRAT_BrGeometryDbManager
#include "BrGeometryDbManager.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#ifndef BRAT_BrCalibration
#include "BrCalibration.h"
#endif
#ifndef BRAT_BrBbRdo
#include "BrBbRdo.h"
#endif
#ifndef BRAT_BrBbVertex
#include "BrBbVertex.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrTofDig
#include "BrTofDig.h"
#endif
#ifndef BRAT_BrTofRdo
#include "BrTofRdo.h"
#endif
#ifndef BRAT_BrDataObject
#include "BrDataObject.h"
#endif
#ifndef BRAT_BrUnits
#include "BrUnits.h"
#endif
#ifndef BRAT_BrSpectrometerTracks
#include "BrSpectrometerTracks.h"
#endif
#ifndef BRAT_BrTofTrackMatch
#include "BrTofTrackMatch.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif

//___________________________________________________________
ClassImp(BrTofRdoModule);

//___________________________________________________________
 BrTofRdoModule::BrTofRdoModule()
  : BrModule()
{
  // Default Constructor. (DO NOT USE)
  fTofParams   = 0;
  fCalibration   = 0;

  SetDefaultParameters();
}

//___________________________________________________________
 BrTofRdoModule::BrTofRdoModule(Char_t *Name,Char_t *Title)
  : BrModule(Name,Title)
{
  fTofParams   = 0;
  fCalibration   = 0;

  if (strcmp(GetName(), "TOF1") == 0) {
    fInFfs = kTRUE;
    fDist  = 900;
  }
  else if (strcasecmp(GetName(), "TOF2") == 0) { 
    fInBfs = kTRUE;
    fDist  = 1800;
  }
  else { 
    fInMrs = kTRUE;
    fDist  = 433.5;
  }

  SetDefaultParameters();
}

//___________________________________________________________
 void BrTofRdoModule::SetDefaultParameters()
{
  // default parameters, called in constructors

  SetMaxTdc();          // default value for hit selection: 4000
  SetEnergyThreshold(); // default value is 0.7
  SetMaxEnergy();       // default is 1.3 
  SetNtuple();          // default is kFALSE
  SetUseSlewing();      // default is kTRUE
  SetUseBbVertex();     // default is kTRUE

  SetUseNewBb();        // default is false        
  fVtxName = "BB";
}

//___________________________________________________________
 BrTofRdoModule::~BrTofRdoModule()
{
  // destructor
}

//___________________________________________________________
 void BrTofRdoModule::DefineHistograms()
{
  // define histograms here
  TDirectory* saveDir = gDirectory;
  
  TDirectory* dir = gDirectory->mkdir(Form("%sRdo",GetName()));
  dir->cd();
  
  // histos
  if (fNtuple)
    fTofTree = 
      new TNtuple(Form("%sRdoTree", GetName()),
		  Form("%s rdo hits tree", GetName()), 
		  "tof:tTime:bTime:bEner:tEner:slat:"
		  "xPos:yPos:tPed:bPed:tAGain:bAGain:"
		  "tTGain:bTGain:timeOff:"
		  "Cs:Time0:Z0:dDelay:trkHit:bbr");
  
  
  fHits = new TH1F("allHits", Form("%s: all hits", GetName()),
		   fTofParams->GetNoSlats(), 0.5, 
		   fTofParams->GetNoSlats() + 0.5);
  
  fSelHits = new TH1F("selectedHits", Form("%s: selected hits", GetName()),
		      fTofParams->GetNoSlats(), 0.5, 
		      fTofParams->GetNoSlats() + 0.5);
  
  
  Float_t tmin = 25;
  if (fInMrs) 
    tmin = 13;
  if (fInBfs) 
    tmin = 55;
  Float_t tmax = tmin+10.;

  fTof = new TH1F("tof", Form("%s: Time of Flight", GetName()), 500, tmin, tmax);


  dir->mkdir("Tof");
  dir->cd("Tof");
  fSTof = new TH1F* [fTofParams->GetNoSlats()];
  for (Int_t i = 0; i < fTofParams->GetNoSlats(); i++)
    fSTof[i] = new TH1F(Form("tof%d",i+1), 
			Form("%s: Time of Flight, slat %d", GetName(), i+1),
			500, tmin, tmax);
  
  dir->cd();

  fSumTof = new TH2F("slatTof", Form("%s: Tof vs Slat", GetName()),
		     fTofParams->GetNoSlats(), 0.5,
		     fTofParams->GetNoSlats() + 0.5, 600, tmin, tmax);

  fHitY = new TH1F("yPos", Form("%s: Hit Y position", GetName()), 300, 
		   -fTofParams->GetSlatHeight()/2. - 5,
		   fTofParams->GetSlatHeight()/2. + 5);

  fHitX = new TH1F("xPos", Form("%s: Hit X position", GetName()),
		   60, -29.5, 30.5);
  fELoss = new TH1F("eLoss", Form("%s: Energy loss", GetName()),
		    300, 0, 4);

  dir->mkdir("ELoss");
  dir->cd("ELoss");
  fSELoss = new TH1F* [fTofParams->GetNoSlats()];
  for (Int_t i = 0; i < fTofParams->GetNoSlats(); i++)
    fSELoss[i] = new TH1F(Form("de%03d", i+1), 
			  Form("%s: Energy loss, slat %d", GetName(),i+1),
			  300, 0, 4);

  dir->cd();
  fSumELoss = new TH2F("eLossSlat", Form("%s: #DeltaE vs Slat", GetName()),
		       fTofParams->GetNoSlats(), 0.5,
		       fTofParams->GetNoSlats() + 0.5, 200, 0, 4);

  dir->mkdir("TvsDE");
  dir->cd("TvsDE");
  fSlew = new TH2F* [fTofParams->GetNoSlats()];
  for (Int_t i = 0; i < fTofParams->GetNoSlats(); i++)
    fSlew[i] = new TH2F(Form("slew%03d", i+1), 
			Form("%s: Slewing effect, slat %d", GetName(),i+1),
			300, 0, 4, 500, tmin, tmax);
  
  gDirectory = saveDir;
  
}

//___________________________________________________________
 void BrTofRdoModule::Init()
{
  //-------------------
  // initialize module
  //-------------------

  if (DebugLevel() > 1)
    cout << "Entering Init of BrTofRdoModule" << endl;
  
  if (fInBfs)
    CheckFsAngles();

  // get important information from managers

  // parameters
  BrParameterDbManager *parDb = BrParameterDbManager::Instance();
  
  if (!parDb) {
    Abort("Init ", "Couldn't instantiate  detector parameter manager");
    return;
  }
  

  fTofParams = 
    (BrDetectorParamsTof*)parDb->
    GetDetectorParameters("BrDetectorParamsTof", (Char_t*)GetName());
  
  // calibration  
  BrCalibrationManager *parMan =
    BrCalibrationManager::Instance();
  
  if (!parMan) {
    Abort("Init", 
	  "Couldn't instantiate calibration manager");
    return;
  }
  
  
  // get the calibration element (created in main)
  fCalibration =
    (BrTofCalibration*)parMan->Register("BrTofCalibration", GetName());
  
  if (!fCalibration) {
    Abort("Init", 
	  "Couldn't get calibration parameters for %s", GetName());
    return;
  }

  // check calibration mode
  BrCalibration::EAccessMode mode = BrCalibration::kRead;
  Int_t nSlats = fTofParams->GetNoSlats();

  // pedestals
  if (!fCalibration->GetAccessMode("topPedestal")      ||
      !fCalibration->GetAccessMode("botPedestal")) {
    fCalibration->Use("topPedestal", mode, nSlats); 
    fCalibration->Use("botPedestal", mode, nSlats); 
  }
  
  // adc gain
  if (!fCalibration->GetAccessMode("topAdcGain")      ||
      !fCalibration->GetAccessMode("botAdcGain")) {
    fCalibration->Use("topAdcGain", mode, nSlats);  
    fCalibration->Use("botAdcGain", mode, nSlats);  
  }
  
  // tdc gain
  if (!fCalibration->GetAccessMode("topTdcGain")      ||
      !fCalibration->GetAccessMode("botTdcGain")) {
    fCalibration->Use("topTdcGain", mode, nSlats);  
    fCalibration->Use("botTdcGain", mode, nSlats);  
  }

  // delta delay
  if (!fCalibration->GetAccessMode("effSpeedOfLight")      ||
      !fCalibration->GetAccessMode("deltaDelay")) {
    fCalibration->Use("effSpeedOfLight", mode, nSlats);
    fCalibration->Use("deltaDelay", mode, nSlats);  
  }

  // time offsets
  if (!fCalibration->GetAccessMode("timeOffset")) {
    fCalibration->Use("timeOffset", mode, nSlats);
  }

  if(fUseSlewing){
    // slewing parameters
    if (!fCalibration->GetAccessMode("topSlewPar1") ||
	!fCalibration->GetAccessMode("topSlewPar2") ||
	!fCalibration->GetAccessMode("botSlewPar1") ||
	!fCalibration->GetAccessMode("botSlewPar2")) {
      fCalibration->Use("topSlewPar1", mode, nSlats);  
      fCalibration->Use("topSlewPar2", mode, nSlats);  
      fCalibration->Use("botSlewPar1", mode, nSlats);  
      fCalibration->Use("botSlewPar2", mode, nSlats);  
    }
  }
  
  // if TOFW, need panel parameters
  if (fInMrs) {
    fPanelPar = new BrDetectorParamsTof * [fTofParams->GetNoPanels()];
    
    for (Int_t p = 0; p < fTofParams->GetNoPanels(); p++)
      fPanelPar[p] = 
	(BrDetectorParamsTof*)parDb->
	GetDetectorParameters("BrDetectorParamsTof", 
			      Form("TFP%d", p+1));
  }
  
  if (fUseNewBb)
    fVtxName = "BB Vertex";
}


//____________________________________________________________________
 void BrTofRdoModule::CheckFsAngles() 
{
  // check if FFS and BFS are aligned
  BrRunInfoManager* runMan = BrRunInfoManager::Instance();
  const BrRunInfo* run = runMan->GetCurrentRun();

  if (run->GetRunNo() == -1) {
    Abort("Begin", "Current run number is -1. Aborting...");
    return;
  } 
  
  if (Verbose()) {
    cout << " ---> FFS at " << run->GetFFSAngle() << " deg. " 
	 << " and BFS at " << run->GetBFSAngle() << " deg. " << endl;
      
    if (run->GetFFSAngle() == run->GetBFSAngle()) {
      fFsLined = kTRUE;
      if (Verbose()) cout << " ---> FS lined up " << endl;
    }
  }
}

//___________________________________________________________
 void BrTofRdoModule::Begin()
{
  // check calibration revision
  if (!fCalibration->RevisionExists("*")) {
    Abort("Begin", "Could not find all need calibration revisions");
    return;
  }
  
  // check FFS and BFS angles
  if (fInBfs)
    CheckFsAngles();

}

//___________________________________________________________
 void BrTofRdoModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // ----------------------------------------------------------
  // Event method to be called once per event.
  // 1- get BB vtx and digtof hits
  // 2- get calibration, set energy, time (ignore slewing for now) 
  //    and position
  // 3- store recontructed hits in event node
  // ----------------------------------------------------------

  SetState(kEvent);
  BrDataObject* bb=0;
  if(fUseBbVertex){
    bb = inNode->GetObject(fVtxName.Data());
    if (!bb) {
      bb = outNode->GetObject(fVtxName.Data());
      if (!bb) {
	if (DebugLevel() > 2)
	  Warning("Event", "No Beam Beam vertex");
	//	return;
      }
    }
  }
  
  Double_t Z0 = 0;
  Double_t T0 = 0;
  Double_t mt = 0;
  Double_t bbr = 0;

  if(fUseBbVertex && bb){
    if (fUseNewBb){
      Z0 = ((BrBbVertex*)bb)->GetZ0();
      T0 = ((BrBbVertex*)bb)->GetTime0();
      mt = ((BrBbVertex*)bb)->GetMethod();
      bbr = ((BrBbVertex*)bb)->GetLeftTime();
    }
    
    else {
      Z0 = ((BrBbRdo*)bb)->GetZ0();
      T0 = ((BrBbRdo*)bb)->GetTime0();
      mt = ((BrBbRdo*)bb)->GetTimeMethod();
    }
    // don't want this vertex (fastest tubes)
    // Beware of dropping events like this
    // comment out for private use
    //    if (mt == 3)
    //  return;
    
    if (Verbose() > 25) {
      if (fUseNewBb)
	((BrBbVertex*)bb)->Print();
      else 
	((BrBbRdo*)bb)->Print();
    }
  }
      

  
  

  //-------------------------------------------------------
  // ------- get match table
  BrDataTable* match = 0;
  if (fInFfs) match = inNode->GetDataTable("FFS Matching");
  if (fInMrs) match = inNode->GetDataTable("MRS Matching");
  if (fInBfs) match = inNode->GetDataTable("BFS Matching");
  if (!match) {
    if (fInFfs) match = outNode->GetDataTable("FFS Matching");
    if (fInMrs) match = outNode->GetDataTable("MRS Matching");
    if (fInBfs) match = outNode->GetDataTable("BFS Matching");
  }
  
  // ------- track table
  BrDataTable* trks  = 0;
  if (fInFfs) trks = inNode->GetDataTable("FfsTracks");
  if (fInMrs) trks = inNode->GetDataTable("MrsTracks");
  if (fInBfs) trks = inNode->GetDataTable("FsTracks");

  if (!trks) {
    if (fInFfs) trks = outNode->GetDataTable("FfsTracks");
    if (fInMrs) trks = outNode->GetDataTable("MrsTracks");
    if (fInBfs) trks = outNode->GetDataTable("FsTracks");
  }
  
  // -------- tof table(s)
  BrDataTable* hits = 
    inNode->GetDataTable(Form("DigTof %s", GetName()));
  if (!hits)
    hits = outNode->GetDataTable(Form("DigTof %s", GetName()));

  // if no hits or tracks
  if (!hits) {
    if (DebugLevel() > 2)
      Warning("Event","No tof table for %s!!", GetName());
    return;
  }

  // prepare output table
  BrTofRdo* recoHits = new BrTofRdo(Form("TofRdo %s", GetName()),
				    Form("TofRdo %s", GetName()));
  
  // ------ tof hits reconstruction   ------
  for(Int_t h = 0; h < hits->GetEntries(); h++) {
    BrTofDig* hit = (BrTofDig*)hits->At(h);
    
    Int_t slat = hit->GetSlatno();
    // ignoring cases where calibration is not good
    if (!fCalibration->ValidCalibration(slat))
      continue;
    
    // ----- get slat calibration
    Float_t    tPed = fCalibration->GetTopPedestal(slat); 
    Float_t    bPed = fCalibration->GetBotPedestal(slat);
    Float_t  tAGain = fCalibration->GetTopAdcGain(slat);
    Float_t  bAGain = fCalibration->GetBotAdcGain(slat);
    Float_t  tTGain = fCalibration->GetTopTdcGain(slat);
    Float_t  bTGain = fCalibration->GetBotTdcGain(slat);
    Float_t    toff = fCalibration->GetTimeOffset(slat);
    Float_t  dDelay = fCalibration->GetDeltaDelay(slat);
    Float_t  cScint = fCalibration->GetEffSpeedOfLight(slat);

    Float_t  slewt1 = 0;
    Float_t  slewt2 = 0;
    Float_t  slewb1 = 0;
    Float_t  slewb2 = 0;

    if (fUseSlewing) {
      slewt1 = fCalibration->GetTopSlewPar1(slat);
      slewt2 = fCalibration->GetTopSlewPar2(slat);
      slewb1 = fCalibration->GetBotSlewPar1(slat);
      slewb2 = fCalibration->GetBotSlewPar2(slat);
      if(slewb2 == BrTofCalibration::kCalException)
	 continue;
    }
    
    // ----- energy loss
    Double_t  tAdc  = hit->GetAdcUp() - tPed ;
    Double_t  tEner = (hit->GetAdcUp() - tPed)/tAGain;
    Double_t  bAdc  = hit->GetAdcDown() - bPed ;
    Double_t  bEner = (hit->GetAdcDown() - bPed)/bAGain;

    // ----- check if energy and tdc ok
    if (tEner < fEnergyThreshold  || bEner < fEnergyThreshold || 
	tEner > fMaxEnergy        || bEner > fMaxEnergy || 
	hit->GetTdcUp() > fMaxTdc || hit->GetTdcDown() > fMaxTdc)
      continue;

    Double_t   ener = 0.5*(tEner + bEner);
    if (HistOn() && tEner > 0 && bEner > 0)
      fHits->Fill(slat);
    
    Double_t  tTime = hit->GetTdcUp()*tTGain;
    Double_t  bTime = hit->GetTdcDown()*bTGain;
    
    // ----- recontruct hit Y pos. 
    Double_t  yPos = (bTime - tTime - dDelay)*cScint/2.;
    
    Int_t    panel = 0;
    Int_t       sl = slat;
    
    BrDetectorParamsTof* tofPar = fTofParams;
    
    // check if TOFW
    if (fInMrs) {
      panel = GetPanelId(slat);
      for (Int_t p = 0; p < panel; p++)
	sl -= fPanelPar[p]->GetNoSlats();
      tofPar = fPanelPar[panel];
    }
    
    // hit X and Z position given by slat pos
    BrVector3D  hitPos = tofPar->GetSlatPos(sl);
    hitPos.SetY(yPos);
    
    // ----- reconstruct time of flight
   
    Double_t tofraw = 0.5 *(tTime+bTime) - toff;
    if(fUseSlewing) {
       Double_t slewCorrection = 0.5 * (slewt2/TMath::Sqrt(tAdc) + slewt1  +
                          slewb2/TMath::Sqrt(bAdc) + slewb1);
       Info(15,"Event","Slewing correction = %f, t1 = %f, t2 = %f",slewCorrection,slewt1,slewt2);
       tofraw -= slewCorrection;
       }
    Double_t tof = tofraw - T0;
    
    Info(14,"Event","Adding hit for slat = %d, tofraw = %f, tof = %f",slat,tof,tofraw);
    recoHits->AddHit(slat, tof, tofraw, ener, hitPos);
    
    // store everything in ntuple (even calibration constants to make
    // some checks in root session)
    if (!HistOn()) 
      continue;
    
    
    fSelHits->Fill(slat);
    fHitX->Fill(hitPos(0));
    fHitY->Fill(hitPos(1));
    fELoss->Fill(0.5*(tEner + bEner));
    fSELoss[slat-1]->Fill(0.5*(tEner + bEner));
    fSumELoss->Fill(slat, 0.5*(tEner + bEner));
    
    // will pick up hits matching tracks to check if slats were aligned
    // in time correctly

    Double_t fixedDistTime = -1;
    Int_t trkHit = 1;

    if (match && trks)
      fixedDistTime = CheckTof(match, trks, slat, tof);

    if (fixedDistTime == -1)
      trkHit = 0;
    
    Debug(10, "Event", " %s: slat %d -> tof at %f = %f",
	  GetName(), slat, fDist, fixedDistTime);
    
    if (trkHit) {
      fTof->Fill(fixedDistTime);
      fSTof[slat-1]->Fill(fixedDistTime);
      fSumTof->Fill(slat, fixedDistTime);
      fSlew[slat-1]->Fill(0.5*(tEner + bEner), fixedDistTime);
    }

    if (fNtuple) {
      const Float_t ntVar[] = { 
	tof, tTime, bTime, tEner, bEner, slat, hitPos(0), hitPos(1),  
	tPed, bPed, tAGain, bAGain, tTGain, bTGain, 
	toff, cScint, T0, Z0, dDelay, trkHit, bbr
      };
      fTofTree->Fill(ntVar);
    }
  }

  // Add to outputNode if more than one
  
  if (recoHits->GetEntries())
    outNode->AddObject(recoHits);
  else
    delete recoHits;
}

//___________________________________________________________
 Double_t BrTofRdoModule::CheckTof(BrDataTable* match, 
				  BrDataTable* trks, 
				  Int_t slat, Double_t tof)
{
  // private method, check tof calibration for hits matching tracks
  
  // --------------------------------------------------------
  // pick up right hit and track
  // we are now sure that the hit will correspond to a valid track
  
  TObject* trk = 0;
  Bool_t hitOk = kFALSE;
  BrTofTrackMatch* pair = 0;
  
  for (Int_t m = 0; m < match->GetEntries(); m++)
    if (((BrTofTrackMatch*)match->At(m))->GetHitId() == slat) {
      hitOk = kTRUE;
      pair = (BrTofTrackMatch*)match->At(m);
      break;
    }
  
  if (!hitOk)
    return -1;
  
  // now select right track proj thanks to the match object
  for (Int_t t = 0; t < trks->GetEntries(); t++) { 
    if (trks->At(t)->GetUniqueID() != pair->GetTrackId())
      continue;
    trk = trks->At(t);
  }
  
  Double_t length;
  if (fInMrs) length = ((BrMrsTrack*)trk)->GetPathLength();
  if (fInFfs) length = ((BrFfsTrack*)trk)->GetPathLength();
  if (fInBfs) length = ((BrFsTrack*)trk)->GetPathLength();
  
  Double_t beta = length/tof;
  return fDist/beta;
}


//___________________________________________________________
 Int_t BrTofRdoModule::GetPanelId(Int_t slat)
{
  // -----------------------------------------
  // returns TOFW panel id given a slat number
  // -----------------------------------------

  Int_t panel = -1;
  Int_t firstSlat[fTofParams->GetNoPanels()];

  // find first slat number of each panel (1 - 21 - 42 - 63 etc)
  firstSlat[0] = 1;
  for (Int_t p = 1; p < fTofParams->GetNoPanels(); p++)
    firstSlat[p] = fPanelPar[p-1]->GetNoSlats() + firstSlat[p-1];
 
  for (Int_t p = 0; p < fTofParams->GetNoPanels(); p++) {
    // last slat of panel p
    Int_t fsl = firstSlat[p] + fPanelPar[p]->GetNoSlats() - 1;
    if (slat >= firstSlat[p] && slat <= fsl) {
      panel = p;
      break;
    }
  }

  return panel;
}
  
//__________________________________________________________________
 void BrTofRdoModule::Print(Option_t* option) const
{
  // Standard information printout. 
  // 
  // Options: See BrModule::Print
  // 
  TString opt(option);
  opt.ToLower();
  
  BrModule::Print(option);
  if (opt.Contains("d"))
    cout << endl
	 << "  Tof Reconstruction Module" << endl
         << "  Original author: Djamel Ouerdane" << endl
         << "  Revisted by:     $Author: hagel $" << endl
         << "  Revision date:   $Date: 2002/08/06 19:07:34 $"   << endl
         << "  Revision number: $Revision: 1.20 $ " << endl
         << endl
         << "*************************************************" << endl;
}

//_______________________________________________________________________________________
//
// $Log: BrTofRdoModule.cxx,v $
// Revision 1.20  2002/08/06 19:07:34  hagel
// Add more verbosity diagnostics
//
// Revision 1.19  2002/06/13 17:26:48  videbaek
// The coding for timing calculation has been modified. The timing offsets are
// always used. The slat/up/down dt part has been removed (thus no dependence on
// effective speed of light in calcualtion of times. The slewings are now
// treated the same in MRS and FS. tested in MRS only though.
// The code also has code to deal with the pp start counters.
//
// Revision 1.18  2002/04/09 02:05:27  ouerdane
// use BrFsTrack instead of BrBfsTrack for H2 rdo check
//
// Revision 1.17  2002/03/20 19:25:02  videbaek
// Upgraded Tof Rdo to save also the raw tof. This is essential for cases where the
// BB counters does not specify the start time and vertex e.g. for FF and MRS
// during the pp run.
//
// Revision 1.16  2001/12/20 15:43:30  ouerdane
// added a poor mans slewing correction for MRS (not optimal but improves PID)
//
// Revision 1.15  2001/12/20 15:37:27  ouerdane
// Modified module to taje into account the slewing correction.
// If the method SetUseSlewing sets the flag to kTRUE,
// The slewing parameters (energy dependence and slewing tail as time offset)
// will be used. If not, only the time offset will be used.
// The default is kFALSE.
// For TOFW, the slewing correction module is not usable since it requires the cherenkov signal
//
// Revision 1.14  2001/12/13 12:45:12  ouerdane
// introduced a poor-man's slewing correction before the proper one is done. Analyses showed it helps!
//
// Revision 1.13  2001/11/19 15:09:40  ouerdane
// updated some time check (histogram level)
//
// Revision 1.12  2001/11/12 18:48:03  ouerdane
// updated in accordance with the new Time offset cal.
//
// Revision 1.11  2001/11/07 10:31:27  ouerdane
// updated rdo module for H2 recontruction and PID
//
// Revision 1.10  2001/11/05 07:06:19  ouerdane
// changed a lot of things to deal with new track classes and BFS Pid
//
// Revision 1.9  2001/10/19 15:47:22  ouerdane
// Removed detector volume object, use now the parameter class method GetSlatPos()
// (correction for the z offsets in H1 and H2 comes naturally with it)
//
// Revision 1.8  2001/10/08 11:29:25  cholm
// Changed to use new DB access classes
//
// Revision 1.7  2001/10/02 20:19:44  ouerdane
// remove member fNoPanels and setter method since it is now part of the
// tof detector parameter class.
//
// Revision 1.6  2001/10/02 01:57:59  ouerdane
// Updated according to the new access mode in BrCalibration
//
// Revision 1.5  2001/08/19 16:00:35  ouerdane
// changed Failure to Abort in Init
//
// Revision 1.4  2001/07/26 13:06:20  ouerdane
// Removed TList private member
// Added diagnostic histograms
// Added option for the ntuple (SetNtuple(Bool_t))
//
// Revision 1.3  2001/06/29 13:52:03  cholm
// Imported TOF classes from Djamel
//
// Revision 1.3  2001/06/21 12:57:17  ouerdane
// Updated to take into account the last modifications of BrTofCalibration
//
// Note: BrTofRdoModule -> added some hardcoded time offset which will
// soon disappear...
//
// Revision 1.2  2001/06/21 12:12:32  ouerdane
// See BrTofCalibration log.
//
// Revision 1.1  2001/06/19 12:46:47  ouerdane
// Added calibration classes and reconstruction classes.
// Brat compiles but these classes haven't been tested in this
// context. Be careful if you use them before I (DO) check if
// all is ok.
//
// Note: some classes are still not included (BrTofSlewingModule,
// BrTofCscintCalModule, BrTofTdcOffsetModule). Will do that after
// Brat2 is available
//

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on by

Validate HTML
Validate CSS