BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//
// $Id: BrTpcCalibration.cxx,v 1.3 2002/08/15 16:15:52 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/08/15 16:15:52 $

//_________________________________________________________________________
//
// BrTpcCalibration is a BRAHMS data class providing storage and 
// access function for TPC specific database parameters.
// This allows modification of the BrCalibrationParams content by the 
// user. 
// Note that all the DB access work is done in the Baseclass. The user
// code tells the Baseclass the name and element of the data
//

#ifndef BRAT_BrTpcCalibration
#include "BrTpcCalibration.h"
#endif

#ifndef BRAT_BrDetectorParamsTPC
#include "BrDetectorParamsTPC.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrDbParameter
#include "BrDbParameter.h"
#endif

//_________________________________________________________________________

ClassImp(BrTpcCalibration);

Char_t* BrTpcCalibration::kDriftVelocity    = "driftVelocity"; 
Char_t* BrTpcCalibration::kInstrumentedRows = "instrumentedRows"; 
Char_t* BrTpcCalibration::kTimeOffset0      = "timeOffset0"; 
Char_t* BrTpcCalibration::kTimeOffset1      = "timeOffset1"; 
Char_t* BrTpcCalibration::kPadStatus        = "padStatus"; 
Char_t* BrTpcCalibration::kPadAdcGain       = "padAdcGain"; 
Char_t* BrTpcCalibration::kPadDeltaTime     = "padDeltaTime"; 

//_________________________________________________________________________
 BrTpcCalibration::BrTpcCalibration()
{
  // Constructor. Set counter and list data members to zero.
  // Don't use this constructor unless you have to and know
  // what you are doing
  // Use BrCalibrationParams(Char_t, Char_t ) instead
  fRows = 0;
  fPadsPrRow = 0;
}

//_________________________________________________________________________
 BrTpcCalibration::BrTpcCalibration(Char_t* name, Char_t* title)
  : BrCalibration(name, title)
{
  // Standard constructor. 
  // the name should be that of the associated detector e.g. "TPM1"
  // Create the Db Objects needed for the specific tables.
  // The names assigned to the variables are only used here and
  // therefore hardcoded

  BrParameterDbManager *gParamDb = BrParameterDbManager::Instance();
  BrDetectorParamsTPC *tpcParams = (BrDetectorParamsTPC*)
    gParamDb->GetDetectorParameters("BrDetectorParamsTPC", GetName());

  fRows = tpcParams->GetNumberOfRows();
  fPadsPrRow = tpcParams->GetPadsprow();
  
  // all
  AddParameterData(kDriftVelocity,    &fDriftVelocity);
  
  //row
  AddParameterData(kInstrumentedRows, &fInstrumentedRows);
  AddParameterData(kTimeOffset0,      &fTimeOffset0);
  AddParameterData(kTimeOffset1,      &fTimeOffset1);
  
  // pad
  AddParameterData(kPadStatus,        &fPadStatus);
  AddParameterData(kPadAdcGain,       &fPadAdcGain);
  AddParameterData(kPadDeltaTime,     &fPadDeltaTime);
}

//_____________________________________________________________________________
 Int_t BrTpcCalibration::GetAccessMode(const Char_t* par) const
{
  // get the accessmode of the parameter data  
  TString s(par);
  Int_t mode = 0;
  
  if (s == kDriftVelocity)
    mode = fDriftVelocity.fAccessMode;
  if (s == kInstrumentedRows)
    mode = fInstrumentedRows.fAccessMode;
  if (s == kTimeOffset0)
    mode = fTimeOffset0.fAccessMode;
  if (s == kTimeOffset1)
    mode = fTimeOffset1.fAccessMode;
  if (s == kPadStatus)
    mode = fPadStatus.fAccessMode;
  if (s == kPadAdcGain)
    mode = fPadAdcGain.fAccessMode;
  if (s == kPadDeltaTime)
    mode = fPadDeltaTime.fAccessMode;
    
  return mode;
}

//_____________________________________________________________________________
Int_t 
 BrTpcCalibration::RowPadTo1D(Int_t row, Int_t pad) const 
{ 
  // Check that row and pad is valid. If valid return the index in the
  // 1 dimensional table
  // rows goes from 1, pads from 0
  if(row < 1 || row > fRows)
    Fatal("RowPadTo1D", "Row %d out of range [1;%d]", row, fRows);
  
  if(pad < 0 || pad >= fPadsPrRow)
    Fatal("RowPadTo1D", "Pad %d out of range [1;%d]", pad, fPadsPrRow);
  
  return (row-1)*fPadsPrRow+pad; 
} 

//_____________________________________________________________________________
Int_t 
 BrTpcCalibration::RowTo1D(Int_t row) const 
{ 
  // Check that row is valid. If valid return row -1 since 
  // rows is 1 based
  if(row < 1 || row > fRows)
    Fatal("RowPadTo1D", "Row %d out of range [1;%d]", row, fRows);
  
  return row-1; 
} 

//_____________________________________________________________________________
Float_t 
 BrTpcCalibration::GetDriftVelocity(void) const {
  if(fDriftVelocity.fAccessMode) 
    return ((Float_t*) fDriftVelocity.fRevision->GetArray())[0] ;
  else
    return 0;
}

//_____________________________________________________________________________
Int_t 
 BrTpcCalibration::GetInstrumentedRows(Int_t row) const {
  if(fInstrumentedRows.fAccessMode) 
    return ((Int_t*) fInstrumentedRows.fRevision->GetArray())[RowTo1D(row)] ;
  else
    return 0;
}

//_____________________________________________________________________________
Short_t 
 BrTpcCalibration::GetPadStatus(Int_t row, Int_t pad) const {
  if(fPadStatus.fAccessMode) 
    return ((Short_t*) fPadStatus.fRevision->GetArray())[RowPadTo1D(row,pad)] ;
  else
    return 0;
}

//_____________________________________________________________________________
Float_t 
 BrTpcCalibration::GetPadAdcGain(Int_t row, Int_t pad) const {
  if(fPadAdcGain.fAccessMode) 
    return ((Float_t*) fPadAdcGain.fRevision->GetArray())[RowPadTo1D(row,pad)] ;
  else
    return 0.0;
}
//_____________________________________________________________________________
Float_t 
 BrTpcCalibration::GetPadDeltaTime(Int_t row, Int_t pad) const {
  if(fPadDeltaTime.fAccessMode) 
    return ((Float_t*) fPadDeltaTime.fRevision->GetArray())[RowPadTo1D(row,pad)] ;
  else
    return 0.0;
}


//_____________________________________________________________________________
Float_t 
 BrTpcCalibration::GetTimeOffset0(Int_t row) const {
  if(fTimeOffset0.fAccessMode) 
    return ((Float_t*) fTimeOffset0.fRevision->GetArray())[RowTo1D(row)] ;
  else
    return 0;
}

//_____________________________________________________________________________
Float_t 
 BrTpcCalibration::GetTimeOffset1(Int_t row) const {
  if(fTimeOffset1.fAccessMode) 
    return ((Float_t*) fTimeOffset1.fRevision->GetArray())[RowTo1D(row)] ;
  else
    return 0;
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetInstrumentedRows(Int_t row, Int_t value) {
  if(fInstrumentedRows.fAccessMode){
    Int_t* array = (Int_t*)fInstrumentedRows.fRevision->GetArray();
    array[RowTo1D(row)]= value;
  }
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetDriftVelocity(Float_t value) {
  if(fDriftVelocity.fAccessMode){
    Float_t* array = (Float_t*)fDriftVelocity.fRevision->GetArray();
    array[0]= value;
  }
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetPadStatus(Int_t row, Int_t pad, Short_t value) {
  if(fPadStatus.fAccessMode){
    Short_t* array = (Short_t*)fPadStatus.fRevision->GetArray();
    array[RowPadTo1D(row,pad)] = value;
  }
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetPadAdcGain(Int_t row, Int_t pad, Float_t value) {
  if(fPadAdcGain.fAccessMode){
    Float_t* array = (Float_t*)fPadAdcGain.fRevision->GetArray();
    array[RowPadTo1D(row,pad)] = value;
  }
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetPadDeltaTime(Int_t row, Int_t pad, Float_t value) {
  if(fPadDeltaTime.fAccessMode){
    Float_t* array = (Float_t*)fPadDeltaTime.fRevision->GetArray();
    array[RowPadTo1D(row,pad)] = value;
  }
}


//_____________________________________________________________________________
void 
 BrTpcCalibration::SetTimeOffset0(Int_t row, Float_t value) {
  if(fTimeOffset0.fAccessMode){
    Float_t* array = (Float_t*)fTimeOffset0.fRevision->GetArray();
    array[RowTo1D(row)]= value;
  }
}

//_____________________________________________________________________________
void 
 BrTpcCalibration::SetTimeOffset1(Int_t row, Float_t value) {
  if(fTimeOffset1.fAccessMode){
    Float_t* array = (Float_t*)fTimeOffset1.fRevision->GetArray();
    array[RowTo1D(row)]= value;
  }
}

//_____________________________________________________________________________
//
// $Log: BrTpcCalibration.cxx,v $
// Revision 1.3  2002/08/15 16:15:52  videbaek
// Commit changes made to be able to apply adcgain correction
// as well as 'timing' offsets. This still has to be worked completely
// out how to calibrate and how to use.
//
// Revision 1.2  2001/11/28 15:18:28  pchristi
// Changed padstatus to short and removed hitwidths.
//
// Revision 1.1  2001/10/12 11:00:47  pchristi
// Added a BrTpcCalibration base class for the TPC calibrations.
//
//

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