BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// Module for tdc delay and cscint calibration: 
//
//  tBot - tTop = yTrack * 2/cscint + Delta_Delay 
//
// fill profiles with tBot - tTop vs yTrack when tracks matching
// hits have a y position between Ymin and Ymax 
// (avoid slat edges where non linearities show up).
// 
// then fit with a 1st degree pol P1:
//  cscint      = inverse slope /2
//  delta delay = P1(y=0)

//____________________________________________________________________
//
// $Id: BrTofDeltaDelayCalModule.cxx,v 1.12 2002/08/15 17:14:38 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/08/15 17:14:38 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrTofDeltaDelayCalModule
#include "BrTofDeltaDelayCalModule.h"
#endif

#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef WIN32
#include <iostream>
#include <fstream>
#include <iomanip>
#else
#include <iostream.h>
#include <fstream.h>
#include <iomanip.h>
#endif

#ifndef ROOT_TH1
#include "TH1.h"
#endif
#ifndef ROOT_TF1
#include "TF1.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrTofDig
#include "BrTofDig.h"
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif
#ifndef BRAT_BrTofTrackMatch
#include "BrTofTrackMatch.h"
#endif
#ifndef BRAT_BrSpectrometerTracks
#include "BrSpectrometerTracks.h"
#endif
#ifndef BRAT_BrCalibration
#include "BrCalibration.h"
#endif

#include <vector>

//____________________________________________________________________
ClassImp(BrTofDeltaDelayCalModule);

//____________________________________________________________________
 BrTofDeltaDelayCalModule::BrTofDeltaDelayCalModule()
  : BrTofCalModule()
{
  // Default constructor. DO NOT USE
  SetState(kSetup);

  SetMinCscint(); // default is 5 
  SetMaxCscint(); // default is 25
}
//____________________________________________________________________
 BrTofDeltaDelayCalModule::BrTofDeltaDelayCalModule(const Char_t* name, 
						   const Char_t* title)
  : BrTofCalModule(name, title)
{
  // Named Constructor
  SetState(kSetup);
  
  SetMinCscint();    // default is 5 
  SetMaxCscint();    // default is 25 
  SetYRange();       // default is -7, 7 cm
}

//____________________________________________________________________
 void BrTofDeltaDelayCalModule::DefineHistograms()
{
  // Define histograms. They are:
  // <fill in here>

  if (GetState() != kInit) {
    Stop("DefineHistograms", "Must be called after Init"); 
    return;  
  }
  
  if (fLoadAscii || fCommitAscii) {
    if (Verbose() > 5)
      Warning("DefineHistograms", "No need to declare histos "
              "since we only load calibration from ascii file"); 
    return;  
  }

  const Int_t nSlats = fParamsTof->GetNoSlats();
  Int_t nPanels=1;
  
  TDirectory* saveDir = gDirectory; 
  fHistDir = gDirectory->mkdir(Form("%sDelays", GetName()));
  fHistDir->cd();
  
  Float_t xmin = -fParamsTof->GetTofWidth()/2.-5;
  Float_t ymin = -fParamsTof->GetSlatHeight()/2.-5;
  Float_t xmax =  fParamsTof->GetTofWidth()/2.+5;
  Float_t ymax =  fParamsTof->GetSlatHeight()/2.+5;
  if (fInMrs) {
    nPanels = fParamsTof->GetNoPanels();
    xmin = -fPanelVol[1]->GetSize()[0]/2.-5;
    xmax =  fPanelVol[1]->GetSize()[0]/2.+5;
  }
  
  //  histograms
  
  fXYTrack = new TH2F* [nPanels];

  for(Int_t k=0; k< nPanels;k++)
    fXYTrack[k] = new TH2F(Form("trkXY_%d",k+1), Form("Y vs X back track projection on tof panel %d",k+1),
			150, xmin, xmax, 150, ymin, ymax);
  
  // histo slat by slat
  fHistDir->mkdir("Slats");
  fHistDir->cd("Slats");
  fDtYTrk = new TProfile* [nSlats];
  for (Int_t i = 1; i <= nSlats; i++) 
    fDtYTrk[i-1] = 
      new TProfile(Form("dty%03d", i), 
		   Form("Slat %d: #DeltaTime vs YTrack", i),
		   100, -fParamsTof->GetSlatHeight()/2.-5, fParamsTof->GetSlatHeight()/2.+5,
		   -10, 10);
  
  fDtYAll = new TH2F("dty", 
		     Form("%s: #DeltaTime vs YTrack all slats", GetName()),
		     200, -fParamsTof->GetSlatHeight()/2.-5, fParamsTof->GetSlatHeight()/2.+5,
		     200, -10, 10);
  fHistDir->cd();
  fSumDelay = new TH1F("delay", 
		       Form("%s: #DeltaDelay (bot - top) in ns", GetName()),
		       nSlats, 0.5, nSlats + 0.5); 
  fSumCscint = new TH1F("csint", 
			Form("%s: Speed of light (cm/ns)", GetName()),
			nSlats, 0.5, nSlats + 0.5); 
  
  gDirectory = saveDir;

}

//____________________________________________________________________
 void BrTofDeltaDelayCalModule::Init()
{
  // Job-level initialisation
  SetState(kInit);
  

  // here is instantialized the calibration parameter and detector
  // parameter objects
  BrTofCalModule::Init();
  
  Int_t nSlats = fParamsTof->GetNoSlats();
  BrCalibration::EAccessMode mode = BrCalibration::kRead;
  
  // check if we want to load cscint and delta delay cal 
  // from ascii file
  if (fLoadAscii) {
    mode = BrCalibration::kTransitional;
    fCalibration->Use("deltaDelay", mode, nSlats);
    fCalibration->Use("effSpeedOfLight", mode, nSlats);
  }
  
  else {
    // check if need calibration already loaded 
    // if not, try to get them from DB
    if (!fCalibration->GetAccessMode("topPedestal")  ||
	!fCalibration->GetAccessMode("botPedestal")) {
      if (DebugLevel() > 3)
	Warning("Init", "Pedestals will be read from the sql DB");
      mode = BrCalibration::kRead;
      fCalibration->Use("topPedestal", mode, nSlats);
      fCalibration->Use("botPedestal", mode, nSlats);
    }
    
    if (!fCalibration->GetAccessMode("topAdcGain")  ||
	!fCalibration->GetAccessMode("botAdcGain")) {
      if (DebugLevel() > 3)
	Warning("Init", "Adc Gains will be read from the sql DB");
      mode = BrCalibration::kRead;
      fCalibration->Use("topAdcGain", mode, nSlats);
      fCalibration->Use("botAdcGain", mode, nSlats);
    }
    
    if (!fCalibration->GetAccessMode("topTdcGain")  ||
	!fCalibration->GetAccessMode("botTdcGain")) {
      if (DebugLevel() > 3)
	Warning("Init", "Tdc Gains will be read from the sql DB");
      mode = BrCalibration::kRead;
      fCalibration->Use("topTdcGain", mode, nSlats);
      fCalibration->Use("botTdcGain", mode, nSlats);
    }
    
    // if we want to save calibration (ascii or DB)
    if (fSaveAscii || fCommitAscii) {
      mode = BrCalibration::kWrite;
      fCalibration->Use("deltaDelay", mode, nSlats);
      fCalibration->Use("effSpeedOfLight", mode, nSlats);
    }
  }
  
  // geometry
  if (!fLoadAscii && !fCommitAscii)
    BrTofCalModule::InitGeo();
}


//____________________________________________________________________
 void BrTofDeltaDelayCalModule::Begin()
{
  // Run-level initialisation
  SetState(kBegin);
  
  if (fLoadAscii) {
    ReadAscii();
    return;
  }
  
  // need histos to do a calibration
  if (!HistOn()) 
    if (!fCommitAscii) {
      Abort("Begin", "Histos are not booked. Cannot proceed"
	    " further.");
      return;
    }
  
  // check if we got parameter revisions:
  if (!fCalibration->RevisionExists("*")) {
    Abort("Begin", "Some calibration revision is missing!");
    return;
  }

  for (Int_t s = 1; s <= fParamsTof->GetNoSlats(); s++) {
    fValidSlat[s - 1] = kTRUE;
    
    Float_t tped = fCalibration->GetTopPedestal(s);
    Float_t bped = fCalibration->GetBotPedestal(s);
    Float_t tag0 = fCalibration->GetTopAdcGain(s);
    Float_t bag0 = fCalibration->GetBotAdcGain(s);
    Float_t ttg  = fCalibration->GetTopTdcGain(s);
    Float_t btg  = fCalibration->GetBotTdcGain(s);
    
    if (!IsValid(tped) || !IsValid(bped) || 
	!IsValid(tag0) || !IsValid(bag0) || 
	!IsValid(ttg) || !IsValid(btg)) {
      fValidSlat[s - 1] = kFALSE;
      if(Verbose()){
	cout << "Not valid Slat from DB " << setw(8) << s << endl;
	if(!IsValid(tped)) cout << "   Top Pedestal" << endl;
	if(!IsValid(bped)) cout << "   Bot Pedestal" << endl;
	if(!IsValid(tag0)) cout << "   Top AdcGain" << endl;
	if(!IsValid(bag0)) cout << "   Bot AdcGain" << endl;
	if(!IsValid(ttg)) cout << "   Top TdcGain" << endl;
	if(!IsValid(btg)) cout << "   Bot TdcGain" << endl;
      }

    }

    fCalibration->SetEffSpeedOfLight(s, BrTofCalibration::kCalException);
    fCalibration->SetDeltaDelay(s, BrTofCalibration::kCalException);
  }
}

//____________________________________________________________________
 void BrTofDeltaDelayCalModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  //---------------------------------------------------------------------//
  //
  //---------------------------------------------------------------------//

  // Per event method
  SetState(kEvent);

  if (fCommitAscii || fLoadAscii)
    return;
  
  // get the the matched tracks and hits table
  BrDataTable* match = 
    inNode->GetDataTable(fMatchName.Data());
  if (!match)
    match = outNode->GetDataTable(fMatchName.Data());

  if (!match) {
    if (Verbose() > 10) Warning("Event", "No tracks matching hits");
    return;
  }

  // get hit table
  BrDataTable* hits = inNode->GetDataTable(Form("DigTof %s", GetName()));
  if (!hits)
    hits = outNode->GetDataTable(Form("DigTof %s", GetName()));
  
  if (!hits) {
    if (Verbose() > 10) Warning("Event", "No global tof hits at all");
    return;
  }

  // get track table
  BrDataTable* tracks = 0;
  if (fInFfs) tracks = inNode->GetDataTable("FfsTracks");
  if (fInBfs) tracks = inNode->GetDataTable("BfsTracks");
  if (fInMrs) tracks = inNode->GetDataTable("MrsTracks");

  if (!tracks) {
    if (fInFfs) tracks = outNode->GetDataTable("FfsTracks");
    if (fInBfs) tracks = outNode->GetDataTable("BfsTracks");
    if (fInMrs) tracks = outNode->GetDataTable("MrsTracks");
  }

  if (!tracks) {
    if (Verbose() > 10) Warning("Event", "No tracks at all");
    return;
  }

  if (Verbose() > 10)
    cout << " -------------------------- " << endl
	 << " Number of tracks         : " << tracks->GetEntries() << endl
         << " Number of tof digits     : " << hits->GetEntries() << endl
	 << " Number of matching pairs : " << match->GetEntries() << endl;
 
  // --------------------------------------------------------

  // now select right hit thanks to the match object
  // we are now sure that the hit will correspond to a valid track  
  
  Int_t nhit = hits->GetEntries();   // number of tof hits in event
  Int_t ntrk = tracks->GetEntries(); // number of tracks in event
  Int_t nmch = match->GetEntries();  // number matching pairs

  BrTofDig*   tofHit[nmch];  // array of good hits
  BrGlbTrack* track [nmch];  // array of tracks

  for (Int_t h = 0; h < nhit; h++) {
    BrTofDig* hit = (BrTofDig*)hits->At(h);
    
    for (Int_t m = 0; m < nmch; m++) {
      BrTofTrackMatch* pair = (BrTofTrackMatch*)match->At(m);
      if (hit->GetSlatno() != pair->GetHitId())
	continue;
      tofHit[m] = hit;
    }
  }
  
  // now select right track proj thanks to the match object

  for (Int_t t = 0; t < ntrk; t++) {
    BrGlbTrack* trk = (BrGlbTrack*)tracks->At(t);
    for (Int_t m = 0; m < nmch; m++) {
      BrTofTrackMatch* pair = (BrTofTrackMatch*)match->At(m);
      if (trk->GetTrackId() != pair->GetTrackId())
	continue;
      track[m] = trk;
    }
  }
  
  //------------------------------------------   
  //    can now proceed with calibration
  //------------------------------------------
  
  for(Int_t h = 0; h < nmch; h++) {
    BrTofTrackMatch* pair = (BrTofTrackMatch*)match->At(h);
    if (!fValidSlat[pair->GetHitId() - 1])
      continue;

    // check calibrated slats
    Int_t   slat = tofHit[h]->GetSlatno();
    Float_t tped = fCalibration->GetTopPedestal(slat);
    Float_t bped = fCalibration->GetBotPedestal(slat);
    Float_t tag0 = fCalibration->GetTopAdcGain(slat);
    Float_t bag0 = fCalibration->GetBotAdcGain(slat);
    Float_t ttg  = fCalibration->GetTopTdcGain(slat);
    Float_t btg  = fCalibration->GetBotTdcGain(slat);

    Double_t tener = (tofHit[h]->GetAdcUp() - tped)/tag0;
    Double_t bener = (tofHit[h]->GetAdcDown() - bped)/bag0;
    
    if (tener < fEnergyThreshold || bener < fEnergyThreshold)
      continue;
    
    Double_t tTime = tofHit[h]->GetTdcUp()*ttg;
    Double_t bTime = tofHit[h]->GetTdcDown()*btg;

    // ------- fill calibration histograms here:
    BrVector3D proj = track[h]->GetProjOnTof(); 
    Int_t panel=0;
    if (fInMrs) {
      panel = pair->GetPanelId();
      fPanelVol[panel]->GlobalToLocal(track[h]->GetProjOnTof(), proj, 0);
    }
    
    fDtYTrk[slat-1]->Fill(proj(1), bTime - tTime);
    fDtYAll->Fill(proj(1), bTime - tTime);
    fXYTrack[panel]->Fill(proj(0), proj(1));
  }
}
 
//____________________________________________________________________
 void BrTofDeltaDelayCalModule::Finish()
{
  // Job-level finalisation
  SetState(kFinish);

  // if load ascii mode
  if (fLoadAscii) 
    return;

  // if commit mode
  if (fCommitAscii) {
    ReadAscii();
    return;
  }
  
  // fit profiles with 1st degree pol. P1
  // inverse slope * 2 = cscint
  // P1(0) = delta delay

  TF1* fit = 0;
  for (Int_t i = 0; i < fParamsTof->GetNoSlats(); i++) {
    Int_t slat = i+1; 
    if (!fValidSlat[slat-1]) {
      fCalibration->SetEffSpeedOfLight(slat, BrTofCalibration::kCalException);
      fCalibration->SetDeltaDelay(slat, BrTofCalibration::kCalException);
      continue;
    }
    
    fDtYTrk[i]->Fit("pol1", "Q0", "", fYmin, fYmax);
    fit = (TF1*)fDtYTrk[i]->GetListOfFunctions()->At(0);
    if (!fit) {
      Warning("Finish", "Something's weird for slat %d", slat);
      fCalibration->SetEffSpeedOfLight(slat, BrTofCalibration::kCalException);
      fCalibration->SetDeltaDelay(slat, BrTofCalibration::kCalException);
      continue;
    }
    
    Double_t cscint = 2/fit->GetParameter(1);
    Double_t dDelay = fit->GetParameter(0);
    fCalibration->SetEffSpeedOfLight(slat, cscint);
    fCalibration->SetDeltaDelay(slat, dDelay);

    if (cscint < fMinCscint || cscint > fMaxCscint) {
      fCalibration->SetEffSpeedOfLight(slat, BrTofCalibration::kCalException);
      fCalibration->SetDeltaDelay(slat, BrTofCalibration::kCalException);
    }

    fSumDelay->Fill(slat, fCalibration->GetDeltaDelay(slat));
    fSumCscint->Fill(slat, fCalibration->GetEffSpeedOfLight(slat));

  }
  
  fCalibration->SetComment("deltaDelay", fComment.Data());
  fCalibration->SetComment("effSpeedOfLight",  fComment.Data());

  // save calibration to ascii file
  if (fSaveAscii)
    SaveAscii();
}

//____________________________________________________________________
 void BrTofDeltaDelayCalModule::SaveAscii() 
{

  // save pedestal to ascii file
  
  BrRunInfoManager* runMan = BrRunInfoManager::Instance();
  const BrRunInfo*     run = runMan->GetCurrentRun();
  if (run->GetRunNo() == -1) {
    Stop("SaveAscii", "RunInfo has run number = -1");
    return;
  }
  
  BrTofCalModule::SaveAscii();
  
  ofstream file(fCalibFile.Data(), ios::out);
  
  file << "****************************************** " << endl;
  file << "*  Calibration for Tof detector " << GetName() << endl;
  file << "*  Delays between bot and top   " << endl;
  file << "*  Effective speed of light     " << endl;
  file << "*  Used events from run " << run->GetRunNo() << endl;
  file << "****************************************** " <<endl;
  file << "*" << endl;    
  file << "* slat  |  delta delay |   cscint  " << endl;
  file << "* ---------------------------------" << endl << endl;
  
  for (Int_t i = 0; i < fParamsTof->GetNoSlats(); i++) {
    Int_t slat = i + 1;
    file << setw(4) << slat
	 << setw(15) << fCalibration->GetDeltaDelay(slat) 
	 << setw(15) << fCalibration->GetEffSpeedOfLight(slat) 
	 << endl;
  }
  
  file << "* ------------------------------------" << endl << endl;
}


//____________________________________________________________________
 void BrTofDeltaDelayCalModule::ReadAscii() 
{
  
  BrTofCalModule::ReadAscii();
  
  ifstream file(fCalibFile.Data(), ios::in);
  
  if (!file) {
    Stop("ReadAscii", "File %s was not found", fCalibFile.Data());
    return;
  }
  
  Float_t delay, cscint;
  Int_t slat;
  Char_t comment[256];
  
  file.getline(comment, 256);
  while(comment[0] == '*') {
    file.getline(comment, 256);
    if (DebugLevel() > 5)
      cout << comment << endl;
  } 

  for (Int_t i = 1; i <= fParamsTof->GetNoSlats(); i++) {
    file >> slat >> delay >> cscint;
    if (DebugLevel() > 5) 
      cout << setw(4) << slat 
	   << setw(12) << delay 
	   << setw(12) << cscint 
	   << endl;
        
    fCalibration->SetDeltaDelay(slat, delay);
    fCalibration->SetEffSpeedOfLight(slat, cscint);
  }
  
  fCalibration->SetComment("deltaDelay",
			   "Generated by BrTofDeltaDelayCalModule: "
			   "fit dTime vs yTrack and took value at "
			   "yTrack = 0");
  fCalibration->SetComment("effSpeedOfLight",
			   "Generated by BrTofDeltaDelayCalModule: "
			   "fit dTime vs yTrack : cscint = 2/slope");
}

//____________________________________________________________________
 void BrTofDeltaDelayCalModule::Print(Option_t* option) const
{
  // Print module information
  // See BrModule::Print for options.
  // In addition this module defines the Option:
  // <fill in here>

  TString opt(option);
  opt.ToLower(); 
  
  BrModule::Print(option); 
  if (opt.Contains("d")) 
   cout << endl 
         << "  Original author: Djamel Ouerdane" << endl
         << "  Last Modifications: " << endl 
         << "    $Author: videbaek $" << endl  
         << "    $Date: 2002/08/15 17:14:38 $"   << endl 
         << "    $Revision: 1.12 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}

//____________________________________________________________________
//
// $Log: BrTofDeltaDelayCalModule.cxx,v $
// Revision 1.12  2002/08/15 17:14:38  videbaek
// Modified diagnostics panel histograms for TOFW
// now one per panel not all panels summed
//
// Revision 1.11  2002/06/13 17:00:42  videbaek
// Added a while back debug output to this otherwise not in use module.
//
// Revision 1.10  2002/03/21 15:04:25  ouerdane
// added fComment member to base class module and method SetComment so that the user can set comments about the calibration at commit time. Removed mean momentum stuff in slewing cal module
//
// Revision 1.9  2001/12/20 15:29:18  ouerdane
// Added method SetYRange(min, max) to set the fitting range
//
// Revision 1.8  2001/11/05 06:59:43  ouerdane
// changed to deal with new track classes and fixed some bugs
//
// Revision 1.7  2001/10/19 15:27:55  ouerdane
// Changed some verbosity level and fixed small things in Init
//
// Revision 1.6  2001/10/08 11:27:57  cholm
// Changed to use new DB access classes
//
// Revision 1.5  2001/10/02 01:53:28  ouerdane
// Added SetSaveAscii, SetLoadAscii, SetCommitAscii and updated the way parameters are tagged for Use
//
// Revision 1.4  2001/08/19 13:43:13  ouerdane
// very minor change (corrected mispelling)
//
// Revision 1.3  2001/08/09 02:16:10  ouerdane
// New version of this class.
// can now do cscint and delta delay in one pass
// see implementation file for more comments
// no need to use BrTofCscintCalModule
//
// Revision 1.2  2001/07/31 09:21:44  ouerdane
// Removed all references to a TList member, replaced
// by histogram members, declare ntuple if fNtuple is true (set
// via SetNtuple)
//
// Revision 1.1  2001/06/29 13:52:19  cholm
// Imported TOF classes from Djamel
//

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