BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// This module is devoted to the TOF detector pedestal calibration 
// It initializes the calibration parameter (BrTofCalibration)
// Fill histograms during the event loops
// fit them with a gaussian during the finish method and put the
// result in the calibration parameter object
//
// depending on the user's desire, this module can save the stuff into
// an ascii file which can be committed to the DB
// (cf. brahms_app/do_app/tof/mm_scripts)
//
//____________________________________________________________________
//
// $Id: BrTofPedCalModule.cxx,v 1.10 2002/03/21 15:04:25 ouerdane Exp $
// $Author: ouerdane $
// $Date: 2002/03/21 15:04:25 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//

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

#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TF1
#include "TF1.h"
#endif
#ifndef BRAT_BrTofPedCalModule
#include "BrTofPedCalModule.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrTofDig
#include "BrTofDig.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif
#ifndef BRAT_BrRunInfo
#include "BrRunInfo.h"
#endif

//________________________________________________________________________

ClassImp(BrTofPedCalModule);

//________________________________________________________________________
 BrTofPedCalModule::BrTofPedCalModule()
  : BrTofCalModule()
{
   //
   // default constructor
   //
  SetState(kSetup);
  fTAdc = 0;
  fBAdc = 0;

}

//________________________________________________________________________
 BrTofPedCalModule::BrTofPedCalModule(Char_t *Name, Char_t *Title) 
  : BrTofCalModule(Name, Title)
{
   //
   // constructor for named object
   //
  fTAdc = 0;
  fBAdc = 0;

  SetState(kSetup);
  SetWidthLimit();  // default is 10 ADC channels
}

//________________________________________________________________________
 BrTofPedCalModule::~BrTofPedCalModule()
{
   //
   // destructor
   //
}

//_______________________________________________________________________
 void BrTofPedCalModule::Init(){
  // ----------------------------------------------------------------
  // Called once per session
  // Get detectorParameters. Inform the parameterElementManeger about
  // database tables to be used (i.e. filled).
  // ----------------------------------------------------------------

  // Job-level initialisation
  SetState(kInit);

  if(DebugLevel() > 1)
    cout << "Entering BrTofPedCalModule::Init() for " 
	 << GetName() << endl;


  // fCalibration and fParamsTof set in base class
  BrTofCalModule::Init();

  Int_t nSlats = fParamsTof->GetNoSlats();

  BrCalibration::EAccessMode mode = BrCalibration::kRead;

  // check if we want to load pedestals from ascii file
  if (fLoadAscii) {
    mode = BrCalibration::kTransitional;
    fCalibration->Use("topPedestal", mode, nSlats);
    fCalibration->Use("topPedestalWidth", mode, nSlats);
    fCalibration->Use("botPedestal", mode, nSlats);
    fCalibration->Use("botPedestalWidth", mode, nSlats);
  }
  
  // check if we want to save them (tof ascii file or DB)
  else if (fSaveAscii || fCommitAscii) {
    mode = BrCalibration::kWrite;
    fCalibration->Use("topPedestal", mode, nSlats);
    fCalibration->Use("topPedestalWidth", mode, nSlats);
    fCalibration->Use("botPedestal", mode, nSlats);
    fCalibration->Use("botPedestalWidth", mode, nSlats);
  }
}

//________________________________________________________________________
 void BrTofPedCalModule::DefineHistograms()
{


  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;  
  }
    
  TDirectory* saveDir = gDirectory;
  fHistDir  = gDirectory->mkdir(Form("%s_Pedestals", GetName())) ;
  
  // ***** histograms section
  
  // create a dir for this detector
  
  fHistDir->cd();

  fHistDir->mkdir("Top");
  fHistDir->mkdir("Bot");
  
  // create histos and store them into list
  
  fHistDir->cd("Top");
  fTAdc = new TH1F* [fParamsTof->GetNoSlats()];

  for(Int_t s = 1; s <= fParamsTof->GetNoSlats(); s++)
    fTAdc[s-1] = new TH1F(Form("top%03d", s),
			  Form("%s: Pedestal top tube %d", GetName(), s),
			  500, 0, 500);
  
  fHistDir->cd("Bot");
  fBAdc = new TH1F* [fParamsTof->GetNoSlats()];

  for(Int_t s = 1; s <= fParamsTof->GetNoSlats(); s++)
    fBAdc[s-1] = new TH1F(Form("bot%03d", s),
			  Form("%s: Pedestal bot tube %d", GetName(), s),
			  500, 0, 500);
  
  
  fHistDir->cd();
  
  // create summaries:
  fTSum = new TH1F("tpeds", Form("%s Top tube Pedestals", GetName()),
		   fParamsTof->GetNoSlats(), 
		   0.5, fParamsTof->GetNoSlats()+.5);
  fTSum->SetMarkerStyle(25);
  fTSum->SetMarkerColor(2);
  fTSum->SetMarkerSize(0.9);
  
  fBSum = new TH1F("bpeds", Form("%s Bot tube Pedestals", GetName()), 
		   fParamsTof->GetNoSlats(), 
		   0.5, fParamsTof->GetNoSlats()+.5);
  fBSum->SetMarkerStyle(25);
  fBSum->SetMarkerColor(4);
  fBSum->SetMarkerSize(0.9);


  fTSumW = new TH1F("twidth", 
		    Form("%s Top tube pedestal widths", GetName()), 
		    fParamsTof->GetNoSlats(), 
		    0.5, fParamsTof->GetNoSlats()+.5);
  fTSumW->SetMarkerStyle(25);
  fTSumW->SetMarkerColor(2);
  fTSumW->SetMarkerSize(0.9);

  
  fBSumW = new TH1F("bwidth", 
		    Form("%s Bot tube pedestal widths", GetName()), 
		    fParamsTof->GetNoSlats(), 
		    0.5, fParamsTof->GetNoSlats()+.5);
  fBSumW->SetMarkerStyle(25);
  fBSumW->SetMarkerColor(2);
  fBSumW->SetMarkerSize(0.9);

  // comes back to previous dir
  gDirectory = saveDir;
}

//________________________________________________________________________
 void BrTofPedCalModule::Begin()
{

  // Run-level initialisation
  SetState(kBegin);

  if (fLoadAscii) {
    ReadAscii();
    return;
  }
  
  // check if histos are booked
  if (!HistBooked())
    if (!fCommitAscii) {
      Abort("Begin", "must book histograms for calibration to work"); 
      return;
    }
  
  for (Int_t s = 1; s <= fParamsTof->GetNoSlats(); s++) {
    fCalibration->SetTopPedestal(s, BrTofCalibration::kCalException);
    fCalibration->SetBotPedestal(s, BrTofCalibration::kCalException);
    fCalibration->SetTopPedestalWidth(s, BrTofCalibration::kCalException);
    fCalibration->SetBotPedestalWidth(s, BrTofCalibration::kCalException);
  }

}


//________________________________________________________________________
 void BrTofPedCalModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  //
  // Fill the histograms that hold the Pedestal spectra 
  // 
  // Required input e.g. : BrDataTable "DigTOF TOF1"
  //
  // 
  
  // Per event method

  
  SetState(kEvent);
  
  if (fLoadAscii || fCommitAscii)
    return;

  if(DebugLevel() > 1)
    cout << "Entering Event() in BrTofPedCalModule for " 
	 << GetName() << endl;
  
  BrDataTable *hits  = 
    inNode->GetDataTable(Form("DigTof %s", GetName()));

  if (!hits) {
    if (Verbose() > 10)
      Warning("Event", "No hit table for %s", GetName());
    return;
  }
  
  Int_t nhit = hits->GetEntries();
  for(Int_t i = 0; i < nhit; i++) {
    BrTofDig* hit = (BrTofDig*)hits->At(i);
    Int_t slat    = hit->GetSlatno();
    
    if (hit->GetTdcUp() < 4000 || hit->GetTdcDown() < 4000) {
      if (DebugLevel() > 20)
	cout << "Hit is not good for pedestal" << endl;
      continue;
    }
    
    if (DebugLevel() > 50)
      cout << GetName() << ": Slat is " << slat 
	   << " -- top adc: "  << hit->GetAdcUp() 
	   << " -- top tdc: "  << hit->GetTdcUp()
	   << " -- bot adc: "  << hit->GetAdcUp() 
	   << " -- bot tdc: "  << hit->GetTdcUp() << endl;
    
    fTAdc[slat-1]->Fill(hit->GetAdcUp());
    fBAdc[slat-1]->Fill(hit->GetAdcDown());
  }
}

//________________________________________________________________________
 void BrTofPedCalModule::Finish()
{
  //
  // called once per job 
  //
  
  // Job-level finalisation
  SetState(kFinish);
  
  // if load ascii mode
  if (fLoadAscii) 
    return;

  // if commit mode
  if (fCommitAscii) {
    ReadAscii();
    return;
  }
  
  cout << GetName() << ": Fitting pedestal with gaussian..."
       << endl;

  if (Verbose())
    cout << " Tube | FitMean | FitSigma | Mean   |  RMS   |"
	 << " DMean  | DSigma " << endl;
  
  TF1* fit = new TF1("fit", "gaus", 0., 1000.);

  for(Int_t i = 1; i <= fParamsTof->GetNoSlats() ; i++) {
    if (fTAdc[i-1]->GetEntries() != 0) {
      // ---- top tubes 
      Float_t mean     = Float_t(fTAdc[i-1]->GetMean());
      Float_t pedWidth = Float_t(fTAdc[i-1]->GetRMS());
      
      fit->SetParameters(fTAdc[i-1]->GetMaximum(), mean, pedWidth);
      
      fTAdc[i-1]->Fit("fit", "Q0", "", 
		      mean - 2*pedWidth, mean + 2*pedWidth);
      
      Double_t meanFit  = fit->GetParameter(1);    
      Double_t sigmaFit = fit->GetParameter(2);
      if (Verbose())
	cout << "top"  << setw(3)  << i
	     << 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;
      }
      
      fTSum->Fill(i, meanFit);
      fTSumW->Fill(i, sigmaFit);
      
      // set parameters for eventual commit into DB
      if (sigmaFit < fWidthLimit) {
	fCalibration->SetTopPedestal(i, meanFit);
	fCalibration->SetTopPedestalWidth(i, sigmaFit);
      }
    }
    
    else 
      if (Verbose())
	cout << " *** " << GetName() << " Tube " <<  i
	     << " : No statistics to perform a fit *** " << endl;

    /////////////////////////////////////////////////////////
    
    // fit bot tube pedestals with gaussian
    if (fBAdc[i-1]->GetEntries() != 0) {
      Float_t mean     = Float_t(fBAdc[i-1]->GetMean());
      Float_t pedWidth = Float_t(fBAdc[i-1]->GetRMS());
      
      fit->SetParameters(fBAdc[i-1]->GetMaximum(), mean, pedWidth);
      fBAdc[i-1]->Fit("fit", "Q0", "", 
		      mean - 2*pedWidth, mean + 2*pedWidth);
      
      Double_t meanFit  = fit->GetParameter(1); // mean
      Double_t sigmaFit = fit->GetParameter(2); // sigma
      
      if (Verbose())
	cout << "bot"  << setw(3)  << i
	     << 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;
      }
      
      fBSum->Fill(i, meanFit);
      fBSumW->Fill(i, sigmaFit);
      
      // set parameters for eventual commit into DB
      if (sigmaFit < fWidthLimit) {
	fCalibration->SetBotPedestal(i, meanFit);
	fCalibration->SetBotPedestalWidth(i, sigmaFit);
      }
      else 
	if (Verbose())
	  cout << " >>>>>> BAD CAL: Width exceeding user's limit <<<<< n";
    }
    else 
      if (Verbose())
	cout << " *** " << GetName() << " Tube " <<  i
	     << " : No statistics to perform a fit *** " << endl;
    
  }
  // set comments
  
  fCalibration->SetComment("topPedestal","Generated by "
			   "BrTofPedCalModule: fit with a gaussian function");
  fCalibration->SetComment("topPedestalWidth","Generated by "
			   "BrTofPedCalModule: fit with a gaussian function");
  fCalibration->SetComment("botPedestal","Generated by "
			   "BrTofPedCalModule: fit with a gaussian function");
  fCalibration->SetComment("botPedestalWidth","Generated by "
			   "BrTofPedCalModule: fit with a gaussian function");

  if (fSaveAscii)
    SaveAscii();

}

//____________________________________________________________________
 void BrTofPedCalModule::SaveAscii() 
{

  // save pedestal to ascii file
    
  BrRunInfoManager* runMan = BrRunInfoManager::Instance();
  Int_t* runs = runMan->GetRunNumbers();
  Int_t  nrun = runMan->GetNumberOfRuns();

  
  BrTofCalModule::SaveAscii();

  ofstream file(fCalibFile.Data(), ios::out);
  
  file << "****************************************** " << endl;
  file << "*  Calibration for Tof detector " << GetName() << endl;
  file << "*  Pedestal calibration " << endl;
  file << "*     Used events from run(s) ";
  for (Int_t r = 0; r < nrun; r++) file << runs[r] << " ";
  file << endl;
  file << "****************************************** " <<endl;
  file << "*" << endl;    
  file << "* slat |  top ped  | top width  |" 
       << "  bot ped  |  bot width  " << endl;
  file << "*  ------------------------------------" 
       << "------------------------" << endl << endl;
  
  for (Int_t i = 0; i < fParamsTof->GetNoSlats(); i++) {
    Int_t slat = i + 1;
    file << setw(5) << slat << setw(13) 
	 << fCalibration->GetTopPedestal(slat) << setw(13)
	 << fCalibration->GetTopPedestalWidth(slat) << setw(13)
	 << fCalibration->GetBotPedestal(slat) << setw(12)
	 << fCalibration->GetBotPedestalWidth(slat) << endl;
  }
  
  file << "* ------------------------------------" << endl << endl;

  cout << GetName() << " pedestals saved to file " 
       << fCalibFile.Data() << endl << endl;
}


//____________________________________________________________________
 void BrTofPedCalModule::ReadAscii() 
{

  // save pedestal to ascii file
  BrTofCalModule::ReadAscii();

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

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

  Float_t tped, twid, bped, bwid;
  Int_t slat;
  Char_t comment[256];
  
  file.getline(comment, 256);
  while(comment[0] == '*')
    file.getline(comment, 256);

  if (DebugLevel() > 1) 
    cout << " ---- Pedestals and widths top and bot ---- " << endl;
   
  for (Int_t i = 1; i <= fParamsTof->GetNoSlats(); i++) {
    file >> slat >> tped >> twid >> bped >> bwid;
    if (DebugLevel() > 1) 
      cout << " slat" << setw(3) << slat 
	   << setw(12) << tped << setw(12) << twid 
	   << setw(12) << bped << setw(12) << bwid << endl;
    
    fCalibration->SetTopPedestal(slat, tped);
    fCalibration->SetTopPedestalWidth(slat, twid);
    fCalibration->SetBotPedestal(slat, bped);
    fCalibration->SetBotPedestalWidth(slat, bwid);
  }
  
  fCalibration->SetComment("topPedestal", fComment.Data());
  fCalibration->SetComment("botPedestal", fComment.Data());
  fCalibration->SetComment("topPedestalWidth", fComment.Data());
  fCalibration->SetComment("botPedestalWidth", fComment.Data());

}

//____________________________________________________________________
 void BrTofPedCalModule::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:04:25 $"   << endl 
         << "    $Revision: 1.10 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}


//------------------------------------------------

// $Log: BrTofPedCalModule.cxx,v $
// 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/11/05 06:59:44  ouerdane
// changed to deal with new track classes and fixed some bugs
//
// Revision 1.8  2001/10/15 23:16:41  ouerdane
// small changes
//
// Revision 1.7  2001/10/08 11:28:04  cholm
// Changed to use new DB access classes
//
// Revision 1.6  2001/10/03 10:10:59  ouerdane
// minor updates (cosmetic)
//
// 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/07/31 09:21:53  ouerdane
// Removed all references to a TList member, replaced
// by histogram members, declare ntuple if fNtuple is true (set
// via SetNtuple)
//
// Revision 1.3  2001/07/13 15:36:26  videbaek
// Add summary histograms for pedestal widths for bot and bottom
//
// Revision 1.2  2001/06/22 17:39:31  cholm
// Changes t odo with data classes renaming.
//
// Revision 1.1.1.1  2001/06/21 14:55:05  hagel
// Initial revision of brat2
//
// Revision 1.2  2001/06/20 09:12:47  ouerdane
// Changed value of BrTofCalibration::kCalException to 999999
// Removed the old exception in BrTofCalModule
//
// Revision 1.1  2001/06/19 12:46:36  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
//
// Revision 1.2  2001/06/01 10:07:51  ouerdane
// This module now derives from BrTofCalModule
//
// Revision 1.1  2001/05/30 17:20:47  ouerdane
// The TOF detector calibration package is released!
// To have all kind of details about it, please consult the README file.
// For trouble-shooting, bugs, comments, send email to ouerdane@nbi.dk
// or brahms-tof-l@bnl.gov
// Enjoy!
//

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