BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// 
// 

//____________________________________________________________________
//
// $Id: BrTpcDriftFibCalModule.cxx,v 1.6 2002/08/19 09:47:29 pchristi Exp $
// $Author: pchristi $
// $Date: 2002/08/19 09:47:29 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrTpcDriftFibCalModule
#include "BrTpcDriftFibCalModule.h"
#endif
#ifndef ROOT_TGraphErrors
#include "TGraphErrors.h"
#endif

#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif

#ifndef ROOT_TF1
#include "TF1.h"
#endif
#ifndef ROOT_TH1
#include "TH1.h"
#endif
#ifndef ROOT_TH2
#include "TH2.h"
#endif
#ifndef BRAT_BrTrack
#include "BrTrack.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif
#ifndef BRAT_BrCalibration
#include "BrCalibration.h"
#endif
#ifndef BRAT_BrException
#include "BrException.h"
#endif
#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif
#ifndef BRAT_BrTableNames
#include "BrTableNames.h"
#endif
#ifndef BRAT_BrTpcFibDig
#include "BrTpcFibDig.h"
#endif

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

//____________________________________________________________________
ClassImp(BrTpcDriftFibCalModule);

//____________________________________________________________________
 BrTpcDriftFibCalModule::BrTpcDriftFibCalModule()
{
  // Default constructor. DO NOT USE
  SetState(kSetup);
  fValidFibres = 0;
  fFibreYPositions = 0;
  fFibreAdcCuts = 0;
  fFibreName = 0;

  SetFitWidth();
  SetFibreZPosition();
}
//____________________________________________________________________
 BrTpcDriftFibCalModule::BrTpcDriftFibCalModule(const Char_t* name, const Char_t* title)
  : BrTpcCalModule(name, title)
{
  // Named Constructor
  SetState(kSetup);
  fValidFibres = new TArrayI(kMaxFibres);
  fFibreYPositions = new TArrayF(kMaxFibres);
  fFibreAdcCuts = new TArrayF(kMaxFibres);
  fFibreName = 0;

  for(Int_t i = 0; i < kMaxFibres; i++) {
    
    fValidFibres->AddAt(0, i);
    fFibreYPositions->AddAt(0, i);
    fFibreAdcCuts->AddAt(0, i);
  }
  
  SetFitWidth();
  SetFibreZPosition();
}

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

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

  TDirectory* saveDir = gDirectory; 
  fHistDir = gDirectory->mkdir(Form("%sDriftFib", fFibreName)); 
  fHistDir->cd();

  hFibres    = new TH1F*[kMaxFibres];
  hFibresAdc = new TH1F*[kMaxFibres];

  // Make histograms here 
  for(Int_t i = 0; i < kMaxFibres; i++) {
    
    hFibres[i] = new TH1F(Form("%s_Fibre%d", fFibreName, i), 
			  Form("%s : Y track position at z = %f for fibre %d above threshold",
			       fFibreName, fFibreZPosition, i),
			  400, -8, 8);

    hFibresAdc[i] = new TH1F(Form("%s_FibreAdc%d", fFibreName, i), 
			  Form("%s : ADC spectrum for fibre %d ",
			       fFibreName, i),
			    128, 0, 2048);
  }

  hFibresXY = new TH2F(Form("%s_FibreXY", fFibreName),
		       Form("%s : X Y projection fibres",fFibreName),
			    40, -20,20, 120, -8, 8 );
  fDriftVelocityGraph = 0;

  gDirectory = saveDir;
}

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

  // base class initialization (register calibration parameters)
  BrTpcCalModule::Init();
  
  BrCalibration::EAccessMode mode = BrCalibration::kRead;
  
  try{
    // check if we want to load adc gain cal from ascii file
    if (fLoadAscii) {
      
      mode = BrCalibration::kTransitional;
    } else {
      
      mode = BrCalibration::kWrite;
    }
    
    fCalibration->Use(BrTpcCalibration::kDriftVelocity, mode, 1);
  }
  catch(BrException* e) {
    cerr << "BrTpcHitWidthCalModule::Init : " << *e << endl;
    e->Execute();
  }
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::Begin()
{
  // Run-level initialisation
  SetState(kBegin);
  
  if (fLoadAscii) {
    ReadAscii();
    SetDetectorParameters();
    return;
  }
  
  if (!HistBooked())
    if (!fCommitAscii) {
      Abort("Begin", "Need histos for calibration");
      return;
    }
  
  // check if we got parameter revisions:
  if (!fCalibration->RevisionExists(BrTpcCalibration::kDriftVelocity)) {
    Abort("Begin", 
	  "Some calibration revision are missing! BrDbUpdateModule ?");
    return;
  }
  
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::SetDetectorParameters()
{
  // Set drift velocity in BrDetectorParamsTPC
  fParamsTpc->SetDriftVelocity(fCalibration->GetDriftVelocity());
}

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

  if (fCommitAscii || fLoadAscii)
    return;
  
  BrDataTable *trackTable = 
    inNode->GetDataTable(Form("%s %s", BRTABLENAMES kTpcTrack, GetName()));

  if(!trackTable) {
    
    if(DebugLevel())
      Warning("Event", "Table was not found : %s %s ", BRTABLENAMES kTpcTrack, GetName());
    return;
  }

  const Int_t nTracks = trackTable->GetEntries();
  
  if(nTracks == 0)
    return;

  BrTpcFibDig* fibObject = (BrTpcFibDig*)
    inNode->GetObject(Form("%s %s", BRTABLENAMES kTpcFibDig, fFibreName));

  if(!fibObject) {
    
    Warning("Event", "Object was not found : %s %s ", BRTABLENAMES kTpcFibDig, fFibreName);
    return;
  }

  for(Int_t i = 0; i < kMaxFibres; i++) {

    if(!fValidFibres->At(i))
      continue;
 
    hFibresAdc[i]->Fill(fibObject->GetAdc(i+1));
    if(fibObject->GetAdc(i+1) < fFibreAdcCuts->At(i))
      continue;

    
    for(Int_t j =0; j < nTracks; j++) {
      
      BrTrack* track = (BrTrack*)trackTable->At(j);
      
      Float_t x, y;
      track->GetXYPositionAtZ(x, y, fFibreZPosition);
      
      hFibres[i]->Fill(y);
      hFibresXY->Fill(x,y);
    }
  }
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::Finish()
{
  // Job-level finalisation In the case of calibrating (fSaveAscii)
  // this is where the fits are done and the drift velocity is
  // extracted
  SetState(kFinish);

  cout << "In finish" << GetName() << endl;

  // if load ascii mode
  if (fLoadAscii) 
    return;
  
  // if commit mode
  if (fCommitAscii) {
    ReadAscii();
    return;
  }

  // Calculate number of fit points
  Int_t nFitPoints = 0;
  
  for(Int_t i = 0; i < kMaxFibres; i++)
    if(fValidFibres->At(i))
      nFitPoints++;

  fDriftVelocityGraph = new TGraphErrors(nFitPoints);
  fHistDir->Append(fDriftVelocityGraph);

  Int_t nPointsAdded = 0;
  
  // Loop over 
  for(Int_t i = 0; i < kMaxFibres; i++) {

    if(hFibres[i]->GetEntries() == 0) {
      delete hFibres[i];
      hFibres[i] = 0;
      continue;
    }
    
    const Int_t maxBin = hFibres[i]->GetMaximumBin();
    const Float_t lowRange = hFibres[i]->GetBinLowEdge(maxBin-fFitWidth);
    const Float_t highRange = hFibres[i]->GetBinLowEdge(maxBin+fFitWidth+1);
    
    hFibres[i]->Fit("gaus", "Q0", "", lowRange, highRange);
    
    TF1* func = hFibres[i]->GetFunction("gaus");

    fDriftVelocityGraph->SetPoint(nPointsAdded, 
				  fFibreYPositions->At(i),
				  func->GetParameter(1));
    fDriftVelocityGraph->SetPointError(nPointsAdded,
				       0,
				       func->GetParError(1));

    nPointsAdded++;
  }
  
  cout << "Finished first fits" << endl;

  if(nPointsAdded!=nFitPoints)
    Error("Finish", "Inconcistent number of fit points");

  fDriftVelocityGraph->Fit("pol1", "Q0", "");

  TF1* func = fDriftVelocityGraph->GetFunction("pol1");

  fCalibration->SetDriftVelocity(fParamsTpc->GetDriftVelocity()/
				 func->GetParameter(1));

  cout << "Saving ASCII file" << endl;

  SaveAscii();
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::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;
  }
  
  BrTpcCalModule::SaveAscii();
  
  ofstream file(fCalibFileName.Data(), ios::out);
  
  file << "****************************************** " << endl;
  file << "*  Calibration for Tpc detector " << GetName() << endl;
  file << "*  Drift velocity calibration         " << endl;
  file << "*  Used events from run " << run->GetRunNo() << endl;
  file << "****************************************** " <<endl;
  file << "*" << endl;    
  file << "* Drift velocity   " << endl;
  file << "* --------------------------------" << endl << endl;
  
  file << setw(15) << fCalibration->GetDriftVelocity() << endl;
  
  file << "* ------------------------------------" << endl << endl;
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::ReadAscii() 
{
  // read calibration from file created by this module
  
  BrTpcCalModule::ReadAscii();
  
  ifstream file(fCalibFileName.Data(), ios::in);
  
  if (!file) {
    Stop("ReadAscii", "File %s was not found", fCalibFileName.Data());
    return;
  }
  
  Char_t comment[256];
  
  file.getline(comment, 256);
  while(comment[0] == '*') {
    file.getline(comment, 256);
    if (DebugLevel() > 5)
      cout << comment << endl;
  } 
  
  Float_t driftVelocity;
  
  file >> driftVelocity;
  
  if (DebugLevel() > 5) 
    cout << setw(15) << driftVelocity << endl;
  
  fCalibration->SetDriftVelocity(driftVelocity);
  
  fCalibration->SetComment(BrTpcCalibration::kDriftVelocity,
			   "Generated by BrTpcDriftFibCalModule: "
			   "fit to fibre peaks");
}

//____________________________________________________________________
 void BrTpcDriftFibCalModule::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: Peter H Christiansen" << endl
         << "  Last Modifications: " << endl 
         << "    $Author: pchristi $" << endl  
         << "    $Date: 2002/08/19 09:47:29 $"   << endl 
         << "    $Revision: 1.6 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}

//____________________________________________________________________
//
// $Log: BrTpcDriftFibCalModule.cxx,v $
// Revision 1.6  2002/08/19 09:47:29  pchristi
// Added old changes to calibration modules.
//
// Revision 1.5  2002/08/15 17:16:50  videbaek
// Modified histograms so fibre peaks are inside histo fro all tpc.
//
// Revision 1.4  2002/03/14 21:45:42  videbaek
// Added some additional histograms for checking of ADC distributions
// and xy position dependence of projections onto fibres.
//
// Revision 1.3  2002/01/03 19:51:06  cholm
// Prepared to use BrTableNames class (or perhaps BrDetectorList) for table names
//
// Revision 1.2  2001/11/13 10:21:35  pchristi
// Changed const to enum to make it more compile friendly
//
// Revision 1.1  2001/11/02 13:38:53  pchristi
// Added new module for calibrating drift velocities using the fibres behind
// T2 and in front of T1.
// Updated old modules with small changes.
//
//

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