BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// Module to determind the centrality of a collision, based on the
// Single Detector Energy (SDE) signal (see also BAN 34).  
//
// It expects to find a BrMultRdo object named Rdo<name>, where <name>
// is the name given to the module, in the input node.  It will write
// a BrMultCent object named Cent<name> to the output node. 
// 

//____________________________________________________________________
//
// $Id: BrMultSdeCentModule.cxx,v 1.6 2001/12/07 19:03:19 cholm Exp $
// $Author: cholm $
// $Date: 2001/12/07 19:03:19 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrMultSdeCentModule
#include "BrMultSdeCentModule.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#ifndef BRAT_BrTileCalibration
#include "BrTileCalibration.h"
#endif
#ifndef BRAT_BrSiCalibration
#include "BrSiCalibration.h"
#endif
#ifndef BRAT_BrTileParameters
#include "BrTileParameters.h"
#endif
#ifndef BRAT_BrSiParameters
#include "BrSiParameters.h"
#endif
#ifdef WIN32
#include <iostream.h>
#include <cfloat.h>
#else
#ifndef __IOSTREAM__
#include <iostream>
#endif
#endif
#ifndef __FSTREAM__
#include <fstream>
#endif
#ifndef __IOMANIP__
#include <iomanip>
#endif
#ifndef __CFLOAT__
#include <cfloat>
#endif
#ifndef BRAT_BrMultRdo
#include "BrMultRdo.h"
#endif
#ifndef BRAT_BrDetectorList
#include "BrDetectorList.h"
#endif
#ifndef BRAT_BrMultCent
#include "BrMultCent.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrSiParameters
#include "BrSiParameters.h"
#endif
#ifndef BRAT_BrTileParameters
#include "BrTileParameters.h"
#endif
#ifndef BRAT_BrTileCalibration
#include "BrTileCalibration.h"
#endif
#ifndef BRAT_BrSiCalibration
#include "BrSiCalibration.h"
#endif
#define BR_CENT_MAX_VERTEX 45.

//____________________________________________________________________
ClassImp(BrMultSdeCentModule);

//____________________________________________________________________
 BrMultSdeCentModule::BrMultSdeCentModule()
{
  // Default constructor. DO NOT USE
  SetState(kSetup);
}
//____________________________________________________________________
 BrMultSdeCentModule::BrMultSdeCentModule(const Char_t* name, 
					 const Char_t* title) 
  : BrModule(name, title), BrMultUtil(name)
{
  // Named Constructor
  // Sets the name to look for. See the class description for more
  // details. 
  SetState(kSetup);
}

//____________________________________________________________________
 void BrMultSdeCentModule::DefineHistograms()
{
  // Define histograms. They are:
  //   TH2D    lowVsHigh      high vs low centrality
  // 

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

  TString hName("Cent"); 
  hName += GetName();
  TDirectory* saveDir = gDirectory; 
  TDirectory* histDir = gDirectory->mkdir(hName.Data()); 
  histDir->cd();

  // Make histograms here 
  //    fLowVsMultHisto  = new TH2F("lowvsMult",  
  //  			      "Multiplicity vs Low cut",  
  //  			      100, 0, 3000, 101, 0, 101);
  //    fLowVsMultHisto->SetXTitle("Multiplicity");
  //    fLowVsMultHisto->SetYTitle("Low cut [%]");
  // 
  //    fHighVsMultHisto = new TH2F("highVsMult", 
  //  			      "Multiplicity vs High cut", 
  //  			      100, 0, 3000, 101, 0, 101);
  //    fHighVsMultHisto->SetXTitle("Multiplicity");
  //    fHighVsMultHisto->SetYTitle("High cut [%]");

  fLowVsHighHisto  = new TH2F("lowVsHigh", "Low vs High Cut", 
			     100, 0, 100, 101, 0, 101);
  fLowVsHighHisto->SetXTitle("Low cut [%]");
  fLowVsHighHisto->SetYTitle("High cut [%]");
  gDirectory = saveDir;
}

//____________________________________________________________________
 void BrMultSdeCentModule::Init()
{
  // Initialise the object. Get the parameters for TMA or SMA, and
  // get a handle on the TMA or SMA calibrations.
  SetState(kInit);
  
  BrParameterDbManager* paramManager = BrParameterDbManager::Instance();
  // Tile Detector parameters 
  TString name = GetName();
  
  if (!name.CompareTo("MultTile"))
    fParameters = (BrSiParameters*)paramManager->
      GetDetectorParameters("BrSiParameters",
			    BrDetectorList::GetDetectorName(kBrSi));
  else if (!name.CompareTo("MultSi"))
    fParameters = (BrTileParameters*)paramManager->
      GetDetectorParameters("BrTileParameters",
			    BrDetectorList::GetDetectorName(kBrTile)); 
  else {
    Failure("Init", "I don't know how to deal with '%s'",
	    name.Data());
    return;
  }
  
  // Test that we got it.
  if(!fParameters) {
    Failure("Init", "no DetectorParameters");
    return;
  }
  
#ifndef BR_MULT_CAL_TMP
  // Calibrations 
  BrCalibrationManager *manager = 
    BrCalibrationManager::Instance();
  
  if (!name.CompareTo("MultTile"))
    fCalibration = static_cast<BrTileCalibration*>
      (manager->Register("BrTileCalibration", 
			 BrDetectorList::GetDetectorName(kBrTile)));
  else if (!name.CompareTo("MultSi"))
    fCalibration = static_cast<BrSiCalibration*>
      (manager->Register("BrSiCalibration", 
			 BrDetectorList::GetDetectorName(kBrSi)));
#else 
  if (!name.CompareTo("MultTile"))
    fCalibration = BrTileTmpCalibration::Instance();
  else if (!name.CompareTo("MultSi"))
    fCalibration   = BrSiTmpCalibration::Instance();
#endif 

  if (!fCalibration) {
    Failure("Init", "didn't get calibration from manager");
    return;
  }
  
  fCalibration->Use("centralityFuncs");
}

    


//____________________________________________________________________
 void BrMultSdeCentModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // Per event method
  // Find the appropiate RD), and then use the calibrations to evalute
  // the centrality. 
  SetState(kEvent);

  // Get the reduced data 
  TString tableName("Rdo");
  tableName += GetName();
  BrMultRdo* rdo = (BrMultRdo*)inNode->GetObject(tableName.Data()); 
  if(!rdo) {
    Stop("Event", "couldn't get %s rdo '%s'", 
	 GetName(),  tableName.Data());
    return;
  }

  // Check if array was hit before doing anything. 
  if (rdo->GetArrayHits() < 1) {
    if (Verbose() > 10) 
      Warning("Event", "no hits in array");
    return;
  }

  // Try to get a vertex, using BrMultUtil::FindVertex
  Double_t vtx = FindVertex(inNode, outNode);

  // Check that we're inside the vertecies for which we have
  // calibrations. 
  if (TMath::Abs(vtx) > BR_CENT_MAX_VERTEX) {    
    if (Verbose() > 8)
      Warning("Event", 
	      "Vertex Z-component %6.2f outside of range +/-%6.2f", 
	      vtx, BR_CENT_MAX_VERTEX);
    return;
  }

  // Get the enery 
  Float_t e = rdo->GetArrayEnergy();
  if (DebugLevel() > 10) 
    cout << "Total energy from " << GetName() << ": " << e << endl;

  // Get the number of functions. 
  Int_t    n   = fCalibration->GetNCentralityFunc();
  if (n < 0) {
    Failure("Event", "Couldn't get any centraility functions");
    return;
  }
  
  // Variables to store progress. 
  Float_t  min   = 0;
  Float_t  max   = 100;
  
  for (Int_t i = 0; i < n; i++) {

    // Calculate current cut value
    Float_t cent = 0;
    Float_t cut  = EvalFunc(i, vtx, cent);

    // Test cut 
    if (e > cut && cent < max) max = cent;
    if (e < cut && cent > min) min = cent;

    // Debug stuff
    if (DebugLevel() > 10)
      cout << "  Cut # " << i << ": " << setw(3) << Int_t(cent) 
	   << "%  E: " << cut << endl;
  }
  // Debug stuff
  if (DebugLevel() > 10) 
    cout << "  E=" << e << "max=" << Int_t(max) 
	 << "% min=" << Int_t(min) << " % " << endl;
    
  BrMultCent* centObj = new BrMultCent(Form("Cent%s", GetName()), 
				       Form("Cent%s", GetName()));
  centObj->SetCent(max);
  centObj->SetMult(min);
  outNode->AddObject(centObj);

  
  // Debug stuff
  if (DebugLevel() > 10) 
    centObj->Print();

  if (HistOn())
    //    fLowVsMultHisto->Fill(e, min);
    //    fHighVsMultHisto->Fill(e, max);
    fLowVsHighHisto->Fill(min, max);
}

//____________________________________________________________________
 Float_t BrMultSdeCentModule::EvalFunc(Int_t i, Float_t x, Float_t& cent) 
{
  // Evaluate the ith function, of order o, at x. Returns FLT_MAX in case
  // of error.   
  Int_t order = fCalibration->GetCentralityFuncOrder();
  if (order < 1) 
    return FLT_MAX;

  cent = fCalibration->GetCentralityFuncCut(i);
  if (cent == FLT_MAX) 
    return FLT_MAX;

  Float_t resu  =  fCalibration->GetCentralityFuncPar(i,0);
  for (Int_t j  = 1; j <= order; j++) 
    resu += fCalibration->GetCentralityFuncPar(i,j) * TMath::Power(x,j);

  if (DebugLevel() > 10)
    cout << "Cent: " << cent << " order: " << order << " x: " << x 
	 << " -> " << resu << endl;

  return resu;
}

//____________________________________________________________________
 void BrMultSdeCentModule::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: Christian Holm Christensen" << endl
         << "  Last Modifications: " << endl 
         << "    $Author: cholm $" << endl  
         << "    $Date: 2001/12/07 19:03:19 $"   << endl 
         << "    $Revision: 1.6 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}

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