BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// Module for pedestal calibration of beam-beam counters 
// 

//____________________________________________________________________
//
// $Id: BrBbPedCalModule.cxx,v 1.7 2002/03/21 15:02:16 ouerdane Exp $
// $Author: ouerdane $
// $Date: 2002/03/21 15:02:16 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrBbPedCalModule
#include "BrBbPedCalModule.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TF1
#include "TF1.h"
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif

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

#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrBbDig
#include "BrBbDig.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif

//____________________________________________________________________
ClassImp(BrBbPedCalModule);

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

  SetMaxWidth(); // default is 10  ADC chan
}

//____________________________________________________________________
 BrBbPedCalModule::BrBbPedCalModule(const Char_t* name, const Char_t* title)
  : BrBbCalModule(name, title)
{
  // Named Constructor
  SetState(kSetup);

  fPedSum      = 0;
  fPedWSum     = 0;

  SetMaxWidth();  // default is 10  ADC chan
}

//____________________________________________________________________
 BrBbPedCalModule::~BrBbPedCalModule()
{
  // destructor
}

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

  if (fCommitAscii || fLoadAscii) {
    if (Verbose() > 10)
      Warning("DefineHistograms", "No need to be here if you "
	      "load or commit an ascii calibration");
    return;
  }

  if (GetState() != kInit) {
    Stop("DefineHistograms", "Must be called after Init"); 
    return;  
  }

  TDirectory* saveDir = gDirectory; 
  TDirectory* histDir = gDirectory->mkdir(Form("%s_Pedestals", GetName())); 
  histDir->cd();

  // Make histograms here 
  histDir->cd();
  
  histDir->mkdir("Pedestal");
  
  histDir->cd("Pedestal");
  fPed = new TH1F* [fNoTubes];
  for(Int_t i = 0; i < fNoTubes; i++)
    fPed[i] = new TH1F(Form("ped%02d", i+1), 
		       Form("Pedestal, tube %d", i+1), 500, 0, 500);
  
  histDir->cd();
  
  // create summary profiles:
  fPedSum  = 
    new TH1F(Form("%sPedestals", GetName()), Form("%s Pedestals ", GetName()),
	     fNoTubes, 0.5, fNoTubes+.5);
  fPedWSum = 
    new TH1F(Form("%sPedestal Width", GetName()), 
	     Form("%s Pedestal Width ", GetName()),
	     fNoTubes, 0.5, fNoTubes+.5);
  
  gDirectory = saveDir;
}

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

  //
  // Called once per session
  // Get detectorParameters. Inform the parameterElementManeger about
  // database tables to be used (i.e. filled).

  if(DebugLevel() > 2)
    cout << "Entering BrBbPedCalModule::Init() for " 
	 << GetName() << endl;

  // initializing base class
  BrBbCalModule::Init();
  BrCalibration::EAccessMode mode = BrCalibration::kRead;

  if(fLoadAscii){
    mode = BrCalibration::kTransitional;
    fCalibration->Use("pedestal", mode, fNoTubes);
    fCalibration->Use("pedestalWidth", mode, fNoTubes);
  }
  else if (fSaveAscii || fCommitAscii) {
    mode = BrCalibration::kWrite;
    fCalibration->Use("pedestal",mode, fNoTubes);
    fCalibration->Use("pedestalWidth", mode, fNoTubes);
  }
}

//____________________________________________________________________
 void BrBbPedCalModule::Begin()
{
  // Run-level initialisation
  SetState(kBegin);
  
  if(fLoadAscii){
    ReadAscii();    
    return;
  }
  
  if (fCommitAscii)
    return;
  
  // check if histos are booked
  if (!HistBooked()) {
    Abort("Begin", "MUST book histograms for calibration to work"); 
    return;
  }

  for (Int_t t = 1; t <= fNoTubes; t++) {
    fCalibration->SetPedestal(t, BrBbCalibration::kCalException);
    fCalibration->SetPedestalWidth(t, BrBbCalibration::kCalException);
  }
}

//____________________________________________________________________
 void BrBbPedCalModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // Per event method
  SetState(kEvent);

  if(DebugLevel() > 1)
    cout << "Entering Event() in BrBbPedCalModule for " 
	 << GetName() << endl;

  if (fCommitAscii || fLoadAscii)
    return;

  BrDataTable *hits = inNode->GetDataTable(fDigName);

  if (hits == 0) {
    if (DebugLevel() > 2)
      Warning("Event", "No hit table for %s", GetName());
    return;
  }
  
  Int_t nhit = hits->GetEntries();
  
  for(Int_t i = 0; i < nhit; i++) {
    BrBbDig* hit = (BrBbDig*)hits->At(i);
    Int_t tube = hit->GetTubeNo();
    Int_t j = tube - 1;
    
    if (hit->GetTdc() < 4000) {
      if (DebugLevel() > 10)
	cout << "Hit is not good for pedestal" << endl;
      continue;
    }
    
    fPed[j]->Fill(hit->GetAdc());
  }
}

//____________________________________________________________________
 void BrBbPedCalModule::End()
{
  // Run-level finalisation
  SetState(kEnd);

}

//____________________________________________________________________
 void BrBbPedCalModule::Finish()
{
  // Job-level finalisation
  SetState(kFinish);
  
  if(fLoadAscii)
    return;

  if (fCommitAscii) {
    ReadAscii();
    return;
  }
  
  cout << " Start fitting procedure for " << GetName()
       << " with simple gaussian..."
       << endl;
  
  TF1* fit = new TF1("fit", "gaus", 0., 1000.);

  if (Verbose())
    cout << " Tube | FitMean | FitSigma | Mean   |  RMS   |"
	 << " DMean  | DSigma " << endl;
  
  for(Int_t i = 0; i < fNoTubes ; i++){
    Int_t tube = i + 1;
    if (fPed[i]->GetEntries() == 0)
      continue;

    // fit tube pedestals with gaussian
    
    Float_t mean     = Float_t(fPed[i]->GetMean());
    Float_t pedWidth = Float_t(fPed[i]->GetRMS());
    
    fit->SetParameters(fPed[i]->GetMaximum(), mean, pedWidth);
    fPed[i]->Fit("fit", "Q0", "", mean - 2*pedWidth, mean + 2*pedWidth);
    
    Double_t meanFit  = fit->GetParameter(1);    
    Double_t sigmaFit = fit->GetParameter(2);
    
    fPedSum->Fill(tube, meanFit);
    fPedWSum->Fill(tube, sigmaFit);
    
    if (Verbose())
      cout << setw(3)  << tube 
	   << setw(10) << setprecision(4) << meanFit 
	   << setw(10) << setprecision(4) << sigmaFit
	   << setw(10) << setprecision(4) << mean 
	   << setw(10) << setprecision(4) << pedWidth
	   << setw(10) << setprecision(4) << meanFit - mean
	   << setw(10) << setprecision(4) << sigmaFit - pedWidth 
	   << endl;
    
    if (TMath::Abs(meanFit - mean) > 20) {
      meanFit  = mean;
      sigmaFit = pedWidth;
      if (Verbose())
	cout << " *** Deviation too large, probably bad fit"
	     << " Take histo MEAN and RMS instead" << endl;
    }
    
    // set parameters for eventual commit into DB
    if (sigmaFit < fMaxWidth) {
      fCalibration->SetPedestal(tube, meanFit);
      fCalibration->SetPedestalWidth(tube, sigmaFit);
    }
    else 
      if (Verbose())
	cout << " >>>>>> BAD CAL: Width exceeding user's limit <<<<< n";
  }

  fCalibration->SetComment("pedestal", fComment.Data());
  fCalibration->SetComment("pedestalWidth", fComment.Data());
  if (fSaveAscii)
    SaveAscii();
  
}

//____________________________________________________________________
 void BrBbPedCalModule::SaveAscii() 
{

  // save pedestal to ascii file

  BrRunInfoManager* runMan = BrRunInfoManager::Instance();
  const BrRunInfo*     run = runMan->GetCurrentRun();
  if (run->GetRunNo() == -1) {
    Abort("SaveAscii", "RunInfo has run number = -1");
    return;
  }
  
  BrBbCalModule::SaveAscii();

  ofstream file(fCalibFile.Data(), ios::out);
  
  file << "****************************************** " << endl;
  file << "*  Calibration for BB counter " << GetName() << endl;
  file << "*     Pedestals " << endl;
  file << "*     Used events from run " << run->GetRunNo() << endl;
  file << "****************************************** " <<endl;
  file << "*" << endl;    
  file << "* tube  |  pedestal  | pedestal width " << endl;
  file << "*  ------------------------------------" << endl << endl;
  
  
  for (Int_t i = 0; i < fNoTubes; i++) {
    Int_t tube = i + 1;
    file << setw(7) << tube << setw(12) 
	 << fCalibration->GetPedestal(tube) << setw(12)
	 << fCalibration->GetPedestalWidth(tube) << endl;
  }
  
  file << "* ------------------------------------" << endl << endl;
}

//____________________________________________________________________
 void BrBbPedCalModule::ReadAscii() 
{

  // save pedestal to ascii file
  BrBbCalModule::SaveAscii();

  ifstream file(fCalibFile.Data(), ios::in);

  if (!file) {
    Abort("ReadAscii", "File %s was not found", fCalibFile.Data());
    return;
  }

  Float_t ped, wid;
  Int_t tube;
  Char_t comment[256];
  
  file.getline(comment, 256);
  while(comment[0] == '*')
    file.getline(comment, 256);
  
  for (Int_t i = 1; i <= fNoTubes; i++) {
    file >> tube >> ped >> wid;
    if (DebugLevel() > 5) {
      cout << setw(4) << tube 
	   << setw(12) << ped 
	   << setw(12) << wid << endl;
    }
        
    fCalibration->SetPedestal(tube, ped);
    fCalibration->SetPedestalWidth(tube, wid);
  }

  fCalibration->SetComment("pedestal","Generated by "
			   "BrBbPedCalModule: fit with a gaussian function");
  fCalibration->SetComment("pedestalWidth","Generated by "
			   "BrBbPedCalModule: fit with a gaussian function");

}


//____________________________________________________________________
 void BrBbPedCalModule::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: ouerdane $" << endl  
         << "    $Date: 2002/03/21 15:02:16 $"   << endl 
         << "    $Revision: 1.7 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}

//____________________________________________________________________
//
// $Log: BrBbPedCalModule.cxx,v $
// Revision 1.7  2002/03/21 15:02:16  ouerdane
// added fComment member to base class module and method SetComment so that the user can set comments about the calibration at commit time
//
// Revision 1.6  2001/10/23 20:47:44  ouerdane
// Updated calibration modules for more efficiency, details are mainly technical
//
// Revision 1.5  2001/10/15 22:46:44  ouerdane
// Fixed some small errors
//
// Revision 1.4  2001/10/15 00:35:07  ouerdane
// Updated all Beam-Beam counter calibration modules (like I did for the TOF
// some weeks ago):
//
// Common changes:
//  added methods Set[Save,Commit,Load]Ascii(Bool_t)
//  removed methods Set[ReadFrom,LoadFrom,SaveTo]File
//
// BrBbGainCalModule:
//  cleaned up a lot the code, added diagnostic histograms showing calibrated
//  ADC after the 1st MIP peak was found.
//  Still at a stage where the 1st MIP peak is the only one to be checked.
//  Later on, will add algorithm to find other peaks.
//
// Added BrBbSlewingCalModule in Makefile.am, etc
// The fit function introduced is dT = a + b/sqrt(dE) + c/dE (where dE = ADC/Gain0)
//
// Revision 1.3  2001/10/08 11:27:19  cholm
// Changed to use new DB access classes
//
// Revision 1.2  2001/09/23 01:51:51  videbaek
// Added prototype Gain calibration module. Updated others to deal with
// Load from File. This may be temporary, but does not break any code as far as
// I can tell.
//
// Revision 1.1  2001/07/23 11:40:01  ouerdane
// Added calibration modules for Beam-Beam counters
//  BrBbCalModule:       base class
//  BrBbPedCalModule:    pedestals
//  BrBbTdcGaincalModule tdc gains
//
// Updated Makefile.am, LinkDef.h, Include.h
//

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