BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//----------------------------------------------------------------//
//                                                                //
//	BrBbCalHitsModule                                         // 
//                                                                // 
//      Module class that takes digitized data (BrBbDig)          //
//      and calibration parameters to reconstruct BB  hits        //
//      create CalHits tables                                     // 
//                                                                // 
//----------------------------------------------------------------//

//
// $Id: BrBbCalHitsModule.cxx,v 1.10 2002/04/23 15:40:11 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/04/23 15:40:11 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> 
//
//----------------------------------------------------------------//
//

#include "BrBbCalHitsModule.h"
#include "iostream.h"
#include "fstream.h"
#include "iomanip.h"
#include "TString.h"
#include "TDirectory.h"
#include "BrParameterDbManager.h"
#include "BrCalibrationManager.h"
#include "BrCalibration.h"
#include "BrBbCalHits.h"
#include "BrDataTable.h"
#include "BrBbDig.h"

ClassImp(BrBbCalHitsModule);

//___________________________________________________________
 BrBbCalHitsModule::BrBbCalHitsModule()
  : BrModule()
{
  // Default Constructor. (DO NOT USE)
  SetDefaultParameters();
}

//___________________________________________________________
 BrBbCalHitsModule::BrBbCalHitsModule(const Char_t* Name, const Char_t* Title)
  : BrModule(Name,Title)
{
  SetDefaultParameters();
}

//___________________________________________________________
 void BrBbCalHitsModule::SetDefaultParameters(){
  fParams   = 0;
  fLCalib   = 0;
  fRCalib   = 0;
  fBblCalHits   = 0;
  fBbrCalHits   = 0;

  SetMaxTdc();          // default value for hit selection: 4000
  SetMinTdc();          // default value for hit selection: 10
  SetEnergyThreshold(); // default value is 0.7
  SetTreeOn();          // default is kFALSE
  SetUseOldCal();       // default is kTRUE
  SetLoadDumpCal();     // default is false 
  SetDumpCalFiles();    // default is '0', '0'
  
  fUsedEvt = 0;
  fAccEvt  = 0;
  fLStat   = 0;
  fRStat   = 0;
  
}


//___________________________________________________________
 BrBbCalHitsModule::~BrBbCalHitsModule()
{
  // destructor
  if (fBblCalHits)
    delete fBblCalHits;
  if (fBbrCalHits)
    delete fBbrCalHits;
}

//___________________________________________________________
 void BrBbCalHitsModule::DefineHistograms()
{
  // define histograms here
  TDirectory* saveDir = gDirectory;
  
  TDirectory* dir = gDirectory->mkdir(Form("%sRdo",GetName()));
  dir->cd();
  
  // histos: in all histos, the tube numbering for the right array is
  // shifted right by 50, so that we can gather everything in same
  // histos

  if (fTreeOn) {
    fBbData = 
      new TTree("bbData", "Bb hit information (left + right)");

    cout << &fBblCalHits << " " << &fBbrCalHits << endl;
    fBbData->Branch("bbLeft",  "BrBbCalHits", &fBblCalHits);
    fBbData->Branch("bbRight", "BrBbCalHits", &fBbrCalHits);
  }
  
  fHits = new TH1F("hits", Form("%s calibrated hits", GetName()), 
		   87, 0.5, 87.5); 
  fTof  = new TH2F("tof", Form("%s: Time of Flight", GetName()), 
		   87, 0.5, 87.5, 2000, 0, 100);
  fEner = new TH2F("ener", Form("%s: Energy loss", GetName()), 
		   87, 0.5, 87.5, 600, 0, 12);
  
  fLMult = new TH1F("lmult", "BB Left Multiplicity",  250, 0.5, 250.5);
  fRMult = new TH1F("rmult", "BB Right Multiplicity", 250, 0.5, 250.5);
  fMult  = new TH2F("mult",  "BB Left mult vs Right mult",
		    250, 0.5, 250.5, 250, 0.5, 250.5);
  fFast  = new TH1F("fast",  "BB Fastest tube", 87, 0.5, 87.5);

  gDirectory = saveDir;
  
}

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

  if(DebugLevel() > 1)
    cout << "Entering Init of BrBbCalHitsModule" << endl;
  
  // get important information from managers
  BrParameterDbManager* parDb = BrParameterDbManager::Instance();
  fParams = (BrDetectorParamsBB*)parDb->GetDetectorParameters("BrDetectorParamsBB", "BB");
  if (!fParams) {
    Abort("Init", "Couldn't get detector parameters");
    return;
  }


  BrCalibrationManager *parMan =
    BrCalibrationManager::Instance();

  if (!parMan) {
    Abort("Init", "Couldn't instantiate calibration manager");
    return;
  }
  
  // get the calibration element (created in main)
  fLCalib = (BrBbCalibration*)parMan->Register("BrBbCalibration", "BBL");
  fRCalib = (BrBbCalibration*)parMan->Register("BrBbCalibration", "BBR");
  
  if (!fLCalib || !fRCalib) {
    Abort("Init", "Couldn't get BB calibration parameters");
    return;
  }

  if (!fLoadDumpCal) {
    // BB Left
    if (!fLCalib->GetAccessMode("pedestal"))
      fLCalib->Use("pedestal", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("pedestalWidth"))
      fLCalib->Use("pedestalWidth", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("adcGain0"))
      fLCalib->Use("adcGain0", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("tdcGain"))
      fLCalib->Use("tdcGain", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("adcGapStart"))
      fLCalib->Use("adcGapStart", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("adcGap"))
      fLCalib->Use("adcGap", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("deltaTdc"))
      fLCalib->Use("deltaTdc", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("slewK"))
      fLCalib->Use("slewK", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("slewDt"))
      fLCalib->Use("slewDt", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fLCalib->GetAccessMode("slewP"))
      fLCalib->Use("slewP", BrCalibration::kRead, fParams->GetNoLeftTubes());
    
    
    // BB right
    
    if (!fRCalib->GetAccessMode("pedestal"))
      fRCalib->Use("pedestal", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("pedestalWidth"))
      fRCalib->Use("pedestalWidth", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("adcGain0"))
      fRCalib->Use("adcGain0", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("tdcGain"))
      fRCalib->Use("tdcGain", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("adcGapStart"))
      fRCalib->Use("adcGapStart", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("adcGap"))
      fRCalib->Use("adcGap", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("deltaTdc"))
      fRCalib->Use("deltaTdc", BrCalibration::kRead, fParams->GetNoLeftTubes());
    if (!fRCalib->GetAccessMode("slewK"))
      fRCalib->Use("slewK", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("slewDt"))
      fRCalib->Use("slewDt", BrCalibration::kRead, fParams->GetNoRightTubes());
    if (!fRCalib->GetAccessMode("slewP"))
      fRCalib->Use("slewP", BrCalibration::kRead, fParams->GetNoRightTubes());
  }

  else {
    if (fBblDumpFile == "0" || fBbrDumpFile == "0") {
      Abort("Init", "You MUST set a file name for the calib. files "
	    "you load. Use SetDumpFiles(<bbl file>, <bbr file>)");
      return;
    }

    // BB Left
    fLCalib->Use("pedestal", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("pedestalWidth", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("adcGain0", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("tdcGain", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("adcGapStart", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("adcGap", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("deltaTdc", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("slewK", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("slewDt", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    fLCalib->Use("slewP", BrCalibration::kTransitional, fParams->GetNoLeftTubes());
    
    // BB right
    fRCalib->Use("pedestal", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("pedestalWidth", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("adcGain0", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("tdcGain", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("adcGapStart", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("adcGap", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("deltaTdc", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("slewK", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("slewDt", BrCalibration::kTransitional, fParams->GetNoRightTubes());
    fRCalib->Use("slewP", BrCalibration::kTransitional, fParams->GetNoRightTubes());
  }
}


//___________________________________________________________
 void BrBbCalHitsModule::Begin()
{
  SetState(kBegin);
  
  // ---------------------------------
  // Check if we have revisions
  // ---------------------------------

  if (!fLCalib->RevisionExists("*")) {
    Abort("Begin", "A revision is missing in Beam-Beam left");
    return;
  }
  
  if (!fRCalib->RevisionExists("*")) {
    Abort("Begin", "A revision is missing in Beam-Beam right");
    return;
  }
  
  if (Verbose() > 15)
    cout << " Found a revision for all BB calibration parameters" 
	 << endl;
  
  if (fLoadDumpCal) { 
    ReadDumpCal(fBblDumpFile.Data(), fParams->GetNoLeftTubes(), fLCalib);
    ReadDumpCal(fBbrDumpFile.Data(), fParams->GetNoRightTubes(), fRCalib);
  }
}

//___________________________________________________________
 void BrBbCalHitsModule::ReadDumpCal(const Char_t* filename, Int_t tubes, 
				    BrBbCalibration* cal)
{
  // private: reads the calib files created from an SQL dump
  ifstream file(filename);

  if (file.fail()) {
    Abort("ReadDumpCal", "File %s was not found", filename);
    return;
  }

  Float_t ped, pedw, ag0, gap, gps, tdg, del, sld, slk, slp;
  Int_t tube;
  Char_t comment[256];
  
  file.getline(comment, 256);
  while(comment[0] == '*') {
    file.getline(comment, 256);
    if (DebugLevel() > 5)
      cout << comment << endl;
  } 
  
  if (DebugLevel())
    cout << " ----- Dump calibration " << filename << ": reading " << endl
	 << "* tube "
	 << "|  Pedestal  |  Ped. Width |  Adc gain0  "    
	 << "|  Adc Gap   |   Gap Start  "
	 << "|  Tdc Gain  |  Delta Tdc  |  Slewing K  "
	 << "| Slewing Dt |  Slewing P "  << endl
	 << "* ---------------------------------------" 
	 << "----------------------------" 
	 << "-----------------------------------------" 
	 << "---------------------------" << endl;  
      
  for (Int_t i = 1; i <= tubes; i++) {
    file >> tube >> ped >> pedw >> ag0 >> gap >> gps >> tdg >> del >> slk >> sld >> slp;
    if (DebugLevel()) 
      cout << setw(5)  << tube  
	   << setw(13) << setprecision(6) << ped 
	   << setw(13) << setprecision(6) << pedw 
	   << setw(13) << setprecision(6) << ag0 
	   << setw(13) << setprecision(6) << gap 
	   << setw(13) << setprecision(6) << gps 
	   << setw(13) << setprecision(6) << tdg 
	   << setw(13) << setprecision(6) << del 
	   << setw(13) << setprecision(6) << slk 
	   << setw(13) << setprecision(6) << sld 
	   << setw(13) << setprecision(6) << slp 
	   << endl;
    
    
    cal->SetPedestal(tube, ped);    cal->SetPedestalWidth(tube, pedw);
    cal->SetAdcGain0(tube, ag0);    cal->SetAdcGap(tube, gap);
    cal->SetAdcGapStart(tube, gps); cal->SetTdcGain(tube, tdg);
    cal->SetDeltaTdc(tube, del);    cal->SetSlewDt(tube, sld);
    cal->SetSlewK(tube, slk);       cal->SetSlewP(tube, slp);
  }
  
  if (DebugLevel()) 
    cout << endl;  

  file.close();
}
    
//___________________________________________________________
 void BrBbCalHitsModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // ----------------------------------------------------------
  // Event method to be called once per event.
  // 1- get raw digits
  // 2- select hit meeting conditions set via the setters
  //    (energy threshold, window on TDC)
  // 3- apply calibration and store things into BrBbRdoNew
  // ----------------------------------------------------------

  if (Verbose() > 11)
    cout << " In BrBbCalHitsModule::Event" << endl;

  SetState(kEvent);

  // BB tables
  BrDataTable *lhits = inNode->GetDataTable("DigBB Left");
  BrDataTable *rhits = inNode->GetDataTable("DigBB Right");
  
  if (!lhits || !rhits) {
    if (Verbose() > 10)
      Warning("Event", "One of the BB tables is missing");
    return;
  }

  // prepare output table
  // make title short - it is afterall stired per event.

  fBblCalHits = new BrBbCalHits("BbCalHits Left", "BB Left hits");
  fBbrCalHits = new BrBbCalHits("BbCalHits Right","BB Right hits");

  fUsedEvt++;

  // add to output node  
  Float_t lMult = Reconstruct(lhits, fLCalib, fBblCalHits);
  Float_t rMult = Reconstruct(rhits, fRCalib, fBbrCalHits);
 
  // some stat
  if (lMult) fLStat++;
  if (rMult) fRStat++;

  // FIXE ME: make it strict for now
  if (lMult && rMult) {
    outNode->AddObject(fBblCalHits);
    outNode->AddObject(fBbrCalHits);
    fAccEvt++;
  }
  else {
    delete fBblCalHits;
    delete fBbrCalHits;
    return;
  }
  
  // diagnostic histos  
  if (!HistOn())
    return;

  fLMult->Fill(fBblCalHits->GetArrayMult());
  fRMult->Fill(fBbrCalHits->GetArrayMult());
  fMult->Fill(fBbrCalHits->GetArrayMult(), fBblCalHits->GetArrayMult());

  // if using tree
  if (fTreeOn && fBblCalHits->GetArrayMult() && fBbrCalHits->GetArrayMult())
    fBbData->Fill();
}

//__________________________________________________________________
 Float_t BrBbCalHitsModule::Reconstruct(BrDataTable* tab,
				      BrBbCalibration* calib,
				      BrBbCalHits* rdo)
{
  //------------------------------------
  // Select hits and apply calibration
  // Fill the rdo object 
  // Return a multiplicity estimation
  //------------------------------------

  // shit right tube numbering for histograms
  Int_t shift = 0;
  if (strcmp(calib->GetName(), "BBR") == 0)
    shift = 50;
  
  // array multiplicity
  Float_t mult = 0;
  Int_t   fast = 0;
  Float_t tfast = 99999;

  for(Int_t h = 0; h < tab->GetEntries(); h++) {
    BrBbDig* digit = (BrBbDig*)tab->At(h);
    Int_t tube     = digit->GetTubeNo();
    
    // check if calibration ok for this tube
    if (!calib->ValidCalibration(tube))
      continue;
    
    Float_t ped = calib->GetPedestal(tube);
    Float_t ag0 = calib->GetAdcGain0(tube);
    Float_t gap = calib->GetAdcGap(tube);
    Float_t gps = calib->GetAdcGapStart(tube);
    Float_t tdg = calib->GetTdcGain(tube);
    Float_t del = calib->GetDeltaTdc(tube);
    Float_t sld = calib->GetSlewDt(tube);
    Float_t slk = calib->GetSlewK(tube);
    Float_t slp = calib->GetSlewP(tube);
    
    // get hit information:
    Double_t adc = digit->GetAdc() - ped;
    Double_t tdc = digit->GetTdc();
    
    // check selection conditions
    if (adc/ag0 < fEnergyThreshold)
      continue;
    
    // adc gap correction
    if (adc > gps)
      adc -= gap;
    
    // time calibration
    if(fUseOldCal) {
      tdc -= (slk/TMath::Sqrt(adc) + sld);
      tdc *= tdg;
    }
    else {
      tdc *= tdg;
      tdc = tdc - del - slk/TMath::Sqrt(adc/ag0) -  slp/(adc/ag0) - sld;
    }
    
    // now select descent time
    if (tdc/tdg < fMinTdc || tdc/tdg > fMaxTdc)
      continue;
    
    // fill rdo object    
    mult += adc/ag0;
    rdo->AddHit(tube, tdc, adc/ag0);
    
    // fastest tube
    if (tdc < tfast) {
      tfast = tdc;
      fast  = tube;
    }
    
    // fill histos if on
    if (!HistOn())
      continue;

    fHits->Fill(tube + shift);
    fEner->Fill(tube + shift, adc/ag0);
    fTof ->Fill(tube + shift, tdc);
  }

  // set multiplicity and fastest tube
  
  rdo->SetArrayMult(mult);
  rdo->SetFastest(fast);

  if (HistOn() && fast != 0)
    fFast->Fill(fast + shift);

  if (Verbose() > 25 && mult)
    rdo->Print("a");

  return mult;
}

//__________________________________________________________________
 void BrBbCalHitsModule::Finish() 
{
  SetState(kFinish);
  
  // prints out some stats:
  if (Verbose() > 2) 
    cout << " ----------------------------"   << endl
	 << "  Beam-Beam RDO statistics:   "  << endl
	 << "  Events used     : " << setw(8) << fUsedEvt << endl
	 << "  Accepted events : " << setw(8) << fAccEvt  << endl
	 << "  Efficiency      : " 
	 << setw(8) << Float_t(fAccEvt)/Float_t(fUsedEvt) << endl
	 << "  Left  CalHits   : " << setw(8) << fLStat   << endl
	 << "  Right CalHits   : " << setw(8) << fRStat   << endl
	 << "  Ratio left/right: " 
	 << setw(8) << Float_t(fLStat)/Float_t(fRStat) << endl
	 << " ----------------------------" << endl;
}

  
//__________________________________________________________________
 void BrBbCalHitsModule::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
	 << "  Bb Reconstruction Module" << endl
         << "  Original author: Djamel Ouerdane" << endl
         << "  Revisted by:     $Author: videbaek $" << endl
         << "  Revision date:   $Date: 2002/04/23 15:40:11 $"   << endl
         << "  Revision number: $Revision: 1.10 $ " << endl
         << endl
         << "*************************************************" << endl;
}


////////////////////////////////////////////////////////////////////////////
//
//   $Log: BrBbCalHitsModule.cxx,v $
//   Revision 1.10  2002/04/23 15:40:11  videbaek
//   Remove ios::in in ifstream (for g++3)
//
//   Revision 1.9  2002/03/20 19:23:24  videbaek
//   ifix small return error, and type in comment
//
//   Revision 1.8  2001/12/13 12:45:50  ouerdane
//   increased verbosity level to 25
//
//   Revision 1.7  2001/12/07 21:23:01  videbaek
//   Shorten title names
//
//   Revision 1.6  2001/11/04 01:56:12  ouerdane
//   Added method to be able to load an ascii calibration dumped from the SQL db. An example is in brat/scripts/calib/bb/BbTestLoadDumpCal.C
//
//   Revision 1.5  2001/10/18 19:12:59  ouerdane
//   Added delta time calibration
//
//   Revision 1.4  2001/10/15 03:57:13  ouerdane
//   changed time recontruction according to new slewing corr. function
//
//   Revision 1.3  2001/10/08 11:29:04  cholm
//   Changed to use new DB access classes
//
//   Revision 1.2  2001/10/02 01:57:59  ouerdane
//   Updated according to the new access mode in BrCalibration
//
//   Revision 1.1  2001/09/23 01:40:59  videbaek
//   Added the module to creatre the BbCalHits data. Implemented based on
//   DOs classes.
//

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