BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// Calibration parameter class for primary vertices offsets
// 

//____________________________________________________________________
//
// $Id: BrVertexCalibration.cxx,v 1.3 2002/05/03 13:57:28 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/05/03 13:57:28 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrVertexCalibration
#include "BrVertexCalibration.h"
#endif
#ifndef BRAT_BrCalibration
#include "BrCalibration.h"
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif

//____________________________________________________________________
ClassImp(BrVertexCalibration);

//____________________________________________________________________
 BrVertexCalibration::BrVertexCalibration()
{
  // Default constructor
}

//_________________________________________________________________________
 BrVertexCalibration::BrVertexCalibration(Char_t *name, Char_t *title) :
  BrCalibration(name, title)
{
  // Standard constructor. 
  // the name should be that of the associated detector e.g. "TOF1"
  // Create the Db Objects needed for the specific tables.
  // Beware difference between constructor and Init. 
  
  AddParameterData("tpcVertexOffset",  &fVertexOffset[0]);
  AddParameterData("bbVertexOffset",   &fVertexOffset[1]);
  AddParameterData("zdcVertexOffset",  &fVertexOffset[2]);
  AddParameterData("inelVertexOffset", &fVertexOffset[3]);
}

//_________________________________________________________________________
 BrVertexCalibration::~BrVertexCalibration()
{
  // Default destructor
  // DO NOT USE
}


//_____________________________________________________________________________
 Int_t BrVertexCalibration::GetAccessMode(const Char_t* detector) const
{
  TString det(detector);
  Int_t mode = 0;
  
  if (det == "TPM1")
    mode = fVertexOffset[kTpc].fAccessMode;
  else if (det == "BB")
    mode = fVertexOffset[kBb].fAccessMode;
  else if (det == "ZDC")
    mode = fVertexOffset[kZdc].fAccessMode;
  else if (det == "Inel")
    mode = fVertexOffset[kInel].fAccessMode;
  else
    Warning("GetAccessMode", "Detector "%s" unknown to me!", detector);
  
  return mode;
}

//_____________________________________________________________________________
 void  BrVertexCalibration::SetVertexOffset(EBrVertexDetector det,
					   EBrSubDetector sub, 
					   Float_t offset, Float_t sigma) 
{
  if (!ValidDetector(det))
    return;
  
  if (fVertexOffset[det].fAccessMode){
    Float_t* array = (Float_t*)fVertexOffset[det].fRevision->GetArray();
    array[sub]   = offset;
    array[sub+1] = sigma;
  }
}

//_____________________________________________________________________________
Float_t 
 BrVertexCalibration::GetVertexOffset(EBrVertexDetector det, EBrSubDetector sub) const 
{
  // get offset for given revision
  if (!ValidDetector(det))
    return 0;
      
  if(fVertexOffset[det].fAccessMode )
    return ((Float_t*) fVertexOffset[det].fRevision->GetArray())[sub] ;
  else
    return 0;
}

//_____________________________________________________________________________
Float_t 
 BrVertexCalibration::GetVertexSigma(EBrVertexDetector det, EBrSubDetector sub) const 
{
  // get offset for given revision
  if (!ValidDetector(det))
    return 0;
  
  if(fVertexOffset[det].fAccessMode )
    return ((Float_t*) fVertexOffset[det].fRevision->GetArray())[sub+1] ;
  else
    return 0;
}

//_____________________________________________________________________________
 Bool_t BrVertexCalibration::ValidCalibration(EBrVertexDetector det, EBrSubDetector sub) 
{
  // check if calibration is not tagged as a cal. exception

  if (!ValidDetector(det))
    return kFALSE;

  if (GetVertexOffset(det, sub) == kCalException)
    return kFALSE;
  
  return kTRUE;
}

//_____________________________________________________________________________
 Bool_t BrVertexCalibration::ValidDetector(EBrVertexDetector det) const
{
  
  if (det < kTpc || det > kInel) {
    Warning("ValidDetector", "Detector bit is wrong. Should only be: kTpc, kBb, kZdc or kInel");
    return kFALSE;
  }
  
  return kTRUE;
}


//____________________________________________________________________
//
// $Log: BrVertexCalibration.cxx,v $
// Revision 1.3  2002/05/03 13:57:28  videbaek
// added code for pp Inel vertex
//
// Revision 1.2  2001/11/05 06:47:24  ouerdane
// added sigma info for vertex analysis
//
// Revision 1.1  2001/11/02 10:37:33  ouerdane
// Added class BrVertexCalibration holding the offsets of the different vertices:
//    BB: big, small and fastest tubes in Z (aligned with TPM1)
//  TPM1: Y offset with beam line
//   ZDC: along Z
//
// The calibration is done so far for the BB counters (cf. modules/calib/vertex)
// Like for other calibrations, the result is stored first to an ascii file
// and then can be committed to the DB after a check).
//
// A new entry in the sql DB ahs been created for that purpose:
//  the "detector" is called VERTEX, its type is BEAM. The parameters are:
//
//   tpcVertexOffset 1 element so far (TPM1)
//    bbVertexOffset 3 elements
//   zdcVertexOffset 1 element but I haven't tried anything with it yet.
//
// (cf also the modules/calib/vertex dir. and scripts/calib/bb/BbVtxOffset.C)
//
//

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