BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
//
// Module to create reduce data objects form detector digits for the
// scintilator tiles in the multiplicity array. 
//

//____________________________________________________________________
// $Id: BrTileRdoModule.cxx,v 1.18 2002/05/08 16:31:01 cholm Exp $
// $Author: cholm $
// $Date: 2002/05/08 16:31:01 $
// $Copyright: 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef BRAT_BrTileRdoModule
#include "BrTileRdoModule.h"
#endif
#ifndef BRAT_BrTileDig
#include "BrTileDig.h"
#endif
#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif
#ifndef BRAT_BrTableNames
#include "BrTableNames.h"
#endif
#ifndef BRAT_BrBbRdo
#include "BrBbRdo.h"
#endif
#ifndef BRAT_BrZdcRdo
#include "BrZdcRdo.h"
#endif
#ifndef BRAT_BrVertex
#include "BrVertex.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrDetectorList
#include "BrDetectorList.h"
#endif
#ifndef BRAT_BrTileParameters
#include "BrTileParameters.h"
#endif
#ifndef BRAT_BrTileCalibration
#include "BrTileCalibration.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#include <iostream.h>
#include <iomanip.h>
#include <limits.h>
#include <float.h>

//____________________________________________________________________
ClassImp (BrTileRdoModule);

//____________________________________________________________________
 BrTileRdoModule::BrTileRdoModule (const Char_t *name,
				  const Char_t *title)
  : BrMultRdoModule(name, title)
{
  // Nothing done here, all done in BrMultRdoModule CTOR
  fSingleAdcLow = 0;
  fSingleAdcHigh = 500;
  fSingleAdcComp = 10;

  fBaseCalAdc = 400; // 2^15 (16 bit dual ADC, one bit used for gap)
  fBaseMult   = 200;   //
  fBaseEnergy = .006;    //

  SetAdcGapLimit();
  SetSaturationChannel();
  SetState(kSetup);
}

//____________________________________________________________________
 void BrTileRdoModule::Init()
{
  //
  SetState(kInit);

  // Detector paramters
  BrParameterDbManager* paramManager = BrParameterDbManager::Instance();
  fParameters = (BrTileParameters*)
    paramManager->GetDetectorParameters("BrTileParameters",
					GetName());
  if(!fParameters)
    Failure("Init", "no DetectortParameters");

  // Calibrations
#ifdef BR_MULT_CAL_TMP
  fCalibration = BrTileTmpCalibration::Instance();
#else
  BrCalibrationManager *manager =
    BrCalibrationManager::Instance();

  fCalibration = static_cast<BrTileCalibration*>
    (manager->Register("BrTileCalibration",
		       BrDetectorList::GetDetectorName(kBrTile)));
#endif
  if (!fCalibration)
    Failure("Init", "didn't get calibration from manager");

  fCalibration->Use("pedestal");
  fCalibration->Use("pedestalWidth");
  fCalibration->Use("adcGap");
  fCalibration->Use("adcGain");
  fCalibration->Use("conversionFunc");
  // These are geometry things proper, but put here for now, until a
  // concensus  of how to do the geometry database.
  fCalibration->Use("ringMap");
  fCalibration->Use("rowMap");
  fCalibration->Use("ringPostion");
  fCalibration->Use("tilt");
}


//____________________________________________________________________
 void BrTileRdoModule::Event(BrEventNode *inNode, BrEventNode *outNode)
{
  //
  //
  //
  SetStatus(kOk);
  SetState(kEvent);
  if (Verbose() > 20)
    cout << "Entering BrTileRdoModule::Event" << endl;

  TString tableName;
  tableName = "Dig";
  tableName += BrDetectorList::GetDetectorName(kBrTile);
  BrTileDig* dig = (BrTileDig*)inNode->GetObject(BRTABLENAMES kDigTiles);
  if(!dig) {
    if (Verbose() > 0) // kShowFailures)
      Warning("Event", "couldn't get Tile digits '%s'", tableName.Data());
    return;
  }

  if (DebugLevel() > 5)
    dig->Print();

  // Try to make a Rdo in output node
  tableName = "Rdo";
  tableName += BrDetectorList::GetDetectorName(kBrTile);
  
  BrMultRdo* rdo = new BrMultRdo(tableName.Data(), "RDO data for Tiles");
  if (!rdo) {
    Abort("Event", "couldn't make a BrTileRdo object!");
    return;
  }
  // Add to output node
  outNode->AddObject(rdo);

  // Some temporary variables for bookkeeping and so on
  Int_t     nRows            = fParameters->GetNoRows();
  Int_t     nRings           = fParameters->GetNoRings();
  Int_t     nTiles           = fParameters->GetNoModules();

  // A loop variable, must be decalered here so that MSVC++/KCC
  // doesn't get confused - poor sods, they should stick to ANSI/ISO
  // standard, alias they don't!
  Int_t     i                = 0;

  // Try to get a vertex.
  fCurrentVtxZ = FindVertex(inNode, outNode);

  // For each ring, set the psuedo-rapidity
  for (i = 0; i < nRings; i++) {
    fCurrentRing = i;
    BrMultRdo::BrRing* ring = rdo->AddRing(fCurrentRing);
    fCurrentRingZ  = fCalibration->GetRingPosition(fCurrentRing);
    fCurrentRingR  = fCalibration->GetRingRadius(fCurrentRing);
    Double_t ringScale=0;
    ring->SetEta(CalibrateAvgEta(ringScale ));
  }

  // Calibrate the individual tile data
  CalibrateIndividual(dig, rdo);

  //__________________________________________________________________
  //
  // End of basic calibrations
  //
  // If we didn't get a vertex, we might as well drop the rest,
  // since it depends on the vertex.
  if (fCurrentVtxZ == FLT_MAX) {
    if (Verbose() > 10)
      Warning("Event", "exiting halfways through because no vertex");
    return;
  }

  //__________________________________________________________________
  //
  // Calibrate the ring energies and multiplicities 
  CalibrateRings(rdo);

  //__________________________________________________________________
  //
  // Calibrate the array energy and multiplicity 
  //
  CalibrateArray(rdo); 

  if (DebugLevel() > 25) 
    rdo->Print("air");  
  else if (DebugLevel() > 10 || Verbose() > 20) 
    rdo->Print("a");  
}

//__________________________________________________________________
 void BrTileRdoModule::CalibrateIndividual(BrTileDig* dig, 
					  BrMultRdo* rdo)
{
  // PRIVATE METHOD:
  // Make the caibrated data for the individual detector elements
  // i.e., each tile. 
  //
  // Please note, that the digit data is currently passed as a
  // TObject, rather than a BrTileDig. This is so that we may read old
  // files with BrDigTiles objects in them. It will soon disappear. 

  // A loop variable, must be decalered here so that MSVC++/KCC
  // doesn't get confused - poor sods, they should stick to ANSI/ISO
  // standard, alias they don't! 
  Int_t i = 0;
  
  // Get some parameters on this detector
  Int_t     nRows            = fParameters->GetNoRows(); 
  Int_t     nRings           = fParameters->GetNoRings(); 
  Int_t     nTiles           = fParameters->GetNoModules();

  // For each tile, calibrate the ADC 
  for (i = 0; i < nTiles; i++) {
    fCurrentSingle   = i;
    fCurrentRing   = fCalibration->GetRingMap(fCurrentSingle) ;         

    // If the ring map says this tile is bad, skip to next
    if (fCurrentRing < 0 || fCurrentRing > nRings) {
      if (DebugLevel() > 10) 
	Warning("Event", "Ring map says tile %d isn't on", 
		fCurrentSingle);
      continue;
    }
    BrMultRdo::BrSingle* single = rdo->AddSingle(fCurrentSingle);

    fCurrentRingZ  = fCalibration->GetRingPosition(fCurrentRing);         
    fCurrentRingR  = fCalibration->GetRingRadius(fCurrentRing);
    BrMultRdo::BrRing* ring = rdo->FindRing(fCurrentRing);
    if (!ring) {
      Failure("CalibrateSingle", "No ring for this single"); 
      return;
    }
    
    // Read the raw ADC value from digit 
    Int_t adc = dig->GetAdc(i + 1);

    Info(25 ,"CalibrateIndividual", "ADC %d = %d", fCurrentSingle,
	 adc); 

    // Fill ADC histogram 
    if (HistOn())
      fSingleAdcHisto->Fill(fCurrentSingle+1, adc);

    // If single is saturated, flag it as such. 
    if (adc > fSaturationChannel) {
      // Maybe we should skip this single here? 
      if (DebugLevel() > 5) 
	Warning("CalibrateSingle", "channel %d is saturated (%d > %d)", 
		fCurrentSingle, adc, fSaturationChannel);
      single->SetSaturated(kTRUE);
    }
    
    //_______________________________________________________________
    // 
    // ADC calibrations 
    //
    // Calculate the calibrated ADC value 
    Bool_t   belowThreshold = kFALSE;
    Double_t singleCalAdc =  CalibrateAdc(adc, belowThreshold);
    // Set the calibrated ADC in output. 
    single->SetCalibratedAdc(Int_t(singleCalAdc));
    // Fill CalAdc histogram 
    if (HistOn())
      fSingleCalAdcHisto->Fill(fCurrentSingle+1, singleCalAdc);
    if(!belowThreshold) {
      // Count total hist;
      rdo->SetArrayHits(rdo->GetArrayHits() + 1);
      // Count hits in ring
      ring->SetHits(ring->GetHits() + 1);
      // Accumelate sum calibrated ADC's 
      rdo->SetArrayCalibratedAdc(Int_t(rdo->GetArrayCalibratedAdc() +
				       singleCalAdc));
      // accumelate sum of calibrated ADCs in ring 
      ring->SetCalibratedAdc(Int_t(ring->GetCalibratedAdc() +
				   singleCalAdc)); 
    }
    Info(25 ,"CalibrateIndividual", "Cal ADC %d = %f", fCurrentSingle,
	 singleCalAdc);
    //_______________________________________________________________
    //
    // End of basic calibrations 
    // 
    // If we didn't get a vertex, we might as well drop the rest,
    // since it depends on the vertex. 
    if (fCurrentVtxZ == FLT_MAX) {
      if (DebugLevel() > 10) 
      	cout << "Skipping Mult for single" << endl;
      continue; 
    }
    
    //_______________________________________________________________
    //
    // Energy calibrations 
    // 
    // Calculate the energy deposited in the single 
    Double_t singleEnergy = 
      CalibrateEnergy(singleCalAdc);
    // Correction factor for vertex position with the factor
    // sin(theta). 
    Double_t distZ = fCurrentVtxZ - fCurrentRingZ; 
    Double_t singleEnergyCor = singleEnergy * fCurrentRingR /
      TMath::Sqrt(distZ * distZ + fCurrentRingR * fCurrentRingR);
    // Set the calibrated Energy in output. 
    single->SetEnergy(singleEnergyCor); 
    // Fill Energy histogram 
    if (HistOn()) 
      fSingleEnergyHisto->Fill(fCurrentSingle+1, singleEnergyCor);
    // accumelate sum of calibrated energy in ring 
    // Please note, that this is not the final energy. We need to do
    // outlier and vertex corrections, but we store the value in the
    // RDO module to save memory 
    ring->SetEnergy(ring->GetEnergy() + singleEnergyCor);
    // Dispite the misleading name, we store the sum square of the
    // energy over the rings, in the member
    // BrSingleRdo::fRingMultiplicity, so that we may save some space. 
    if(belowThreshold) 
      continue;
    ring->SetMultiplicity(ring->GetMultiplicity() + 
			  singleEnergyCor * singleEnergyCor);
    //_______________________________________________________________
    //
    // Multiplicity calibrations 
    // 
    Double_t ringEta    = ring->GetEta();
    Double_t singleMult = TMath::Max(0., singleEnergy *
				     CalibrateMultiplicity(ringEta));
    single->SetMultiplicity(singleMult);
    if (HistOn()) 
      fSingleMultHisto->Fill(fCurrentSingle+1, singleMult);
  }

  // Fill histograms 
  if (HistOn()) {
    fArrayCalAdcHisto->Fill(rdo->GetArrayCalibratedAdc());
    fArrayHitsHisto->Fill(rdo->GetArrayHits());
    fVtxVsCalAdcHisto->Fill(fCurrentVtxZ, rdo->GetArrayCalibratedAdc());

    // Fill out hits and calibrated ADCs
    for (i = 0; i < nRings; i++) {
      fRingCalAdcHisto->Fill(i+1, rdo->GetRingCalibratedAdc(i));
      fRingHitsHisto->Fill(i+1, rdo->GetRingHits(i));
    }
    // Reset counter 
  }
}

//__________________________________________________________________
 Double_t BrTileRdoModule::CalibrateArrayMultiplicity()
{
  // PRIVATE METHOD:
  // Correct the sum multiplicity for the vertex position. 
  // This is based in on fits from GEANT, as far as I gather. 
  if (fCurrentVtxZ == FLT_MAX) {
    if (DebugLevel() > 25)
      Warning("CalibrateArrayMultiplicity","A null vertex");
    return 0;
  }
  
  // Get the order of the correcction function
  Int_t order = fCalibration->GetCorrectionFuncOrder();
  if (order < 0) {
    if (DebugLevel() > 25)
      Warning("CalibrateArrayMultiplicity","A null function");
    return 0;
  }

  // Get the valid vertex Z range for the correction function
  Float_t cut = fCalibration->GetCorrectionFuncCut();
  if (cut < 0) {
    if (DebugLevel() > 25)
      Warning("CalibrateArrayMultiplicity", "a null range");
    return 0;
  }
  
  // If we are outside of the valid range in Z, return 0
  if (cut < TMath::Abs(fCurrentVtxZ)) {
    if (DebugLevel() > 25)
      Warning("CalibrateArrayMultiplicity", 
	      "outside range: abs(%f) > %f", 
	      fCurrentVtxZ, cut);
    return 0;
  }

  // Calculate the correction 
  Double_t corr = fCalibration->GetCorrectionFuncPar(0);
  for (Int_t i = 1; i < order + 1; i++) 
    corr += fCalibration->GetCorrectionFuncPar(i) 
      * TMath::Power(fCurrentVtxZ, i); 

  // Show the calibration 
  if (DebugLevel() > 25) {
    cout << GetName() << " correction: (" << fCurrentVtxZ << ") " 
	 << corr << " = " << fCalibration->GetCorrectionFuncPar(0) 
	 << flush; 
    for (Int_t i = 1; i < order + 1; i++) 
      cout << " + " << fCalibration->GetCorrectionFuncPar(i) 
	   << " * " << TMath::Power(fCurrentVtxZ, i) << flush;
    cout << endl;
  }
      
  // if the correction is very small, return 0;
  if (TMath::Abs(corr) < 1e-6)
    return 0;
  
  return 1. / corr;
}

  
//____________________________________________________________________
 void BrTileRdoModule::Print(Option_t* option) const
{ 
  // Module Information method
  //
  // Options: (see also BrModule::Print)
  //  
  TString opt(option);
  opt.ToLower();
  
  BrMultRdoModule::Print(option);
  if (opt.Contains("d")) {    
    cout << endl
         << "  Original author: Christian Holm Christensen" << endl
         << "    $Author: cholm $" << endl 
         << "    $Date: 2002/05/08 16:31:01 $"   << endl
         << "    $Revision: 1.18 $ " << endl 
         << endl
         << "-------------------------------------------------" << endl;
  }
} 


//____________________________________________________________________
//
// $Log: BrTileRdoModule.cxx,v $
// Revision 1.18  2002/05/08 16:31:01  cholm
// Added 3 histograms, one of the sum of the multiplicity in each eta bin,
// one with the number of fills, and a third, the ratio of these two.  Unless
// something tricky is going on, this should be our dN/deta distribution for
// the TMA and SMA.  Also fixed a few warning (probably GCC3 errors) in the
// Si and Tile modules.
//
// Revision 1.17  2002/03/05 20:35:53  sanders
// Fixed problem which kept low hit event being counted in dNdEta
// calculations.
//
// Revision 1.16  2002/01/27 16:13:08  cholm
// Use pass-by-reference ('&') rather than pointers ('*') in the method
// BrMultRdoModule::CalibrateAdc, since that's the proper C++ way to go
// about these things.
//
// Revision 1.15  2002/01/26 16:17:33  sanders
// Modified point where the energy threshold is applied.  Delaying the
// threshold cut until after the energy calibration is necessary for
// checking the energy calibration.
// Also corrected an indexing error that was cutting off the first adc
// channel in the histograms.
//
// Revision 1.14  2002/01/04 21:30:11  sanders
// Added "Setters" for range of raw adc values.
//
// Revision 1.13  2002/01/03 19:52:30  cholm
// Prepared to use BrTableNames class (or perhaps BrDetectorList) for table names
//
// Revision 1.12  2001/12/14 15:41:03  cholm
// Updated for new BrModule::Info method
//
// Revision 1.11  2001/12/07 18:13:33  cholm
// Set better ranges for the ADC, Energy and Multiplicity histograms.
//
// Revision 1.10  2001/11/13 15:51:35  cholm
// REmoved initialisation of calibrations.
//
// Revision 1.9  2001/11/12 14:55:46  sanders
// Added modifications for 200 GeV multiplicity/centrality calibrations.
// Added cfs method to BrZdcRdoModule since this is needed to obtain
// reliable position for low amplitude signals (avoids strong slewing
// correction)
//
// Revision 1.8  2001/10/29 21:01:26  cholm
// Corrected a very bad mistake in the previous commit:  A manager was
// initialised in the module - BAD.  Managers may not be initialised in
// modules, only thier services may be used.
//
// Revision 1.7  2001/10/27 22:58:59  sanders
// Modified files to allow for yr 2 calibration.  The run number is
// now being found using BrRunInfoManager
//
// Revision 1.6  2001/10/08 11:29:19  cholm
// Changed to use new DB access classes
//
// Revision 1.5  2001/09/20 13:44:22  cholm
// Removed the use of temp ASCII files, since not needed anymore
//
// Revision 1.4  2001/08/03 21:36:31  zdc
// The methods to read low-gain ZDC ADC have been added to BrZdcRdoModule
//
// Revision 1.3  2001/07/20 15:18:42  cholm
// A few debug lines, that is it.
//
// Revision 1.2  2001/06/22 17:41:56  cholm
// Change names so that every data class has the same format, so that
// we will not have to worry about that later on. The affected classes
// are:
//
//         BrDigBB             ->        BrBbDig
//         BrRdoBB             ->        BrBbRdo
//         BrDigZDC            ->        BrZdcDig
//         BrRdoZDC            ->        BrZdcRdo
//         BrDigRHIC           ->        BrRichDig
//         BrDigDC             ->        BrDcDig
//         BrDigC1             ->        BrDcC1
// 	BrDigHit            ->	      BrHitDig
// 	BrDigTof	    ->	      BrTofDig
// 	BrTPCSequence	    ->	      BrTpcSequence
// 	BrTPCCluster	    ->	      BrTpcCluster
// 	BrTPCClusterTable   ->	      BrTpcClusterTable
//
// These changes has ofcourse been propegated to the modules as well,
// giving the changes
//
//         BrRdoModuleBB	    ->	      BrBbRdoModule
// 	BrRdoModuleZDC	    ->	      BrZdcRdoModule
// 	BrTPCClusterFinder  ->	      BrTpcClusterFinder
// 	BrTPCSequenceAdder  ->	      BrTpcSequenceAdder
//
// Revision 1.1.1.1  2001/06/21 14:55:09  hagel
// Initial revision of brat2
//
// Revision 1.23  2001/06/01 15:49:36  cholm
// Fix of some default values from -1 to zero (gave bad results)
// BrSiCalibration had a bug in returning the correction function.
// Tile|Rdo modules had some minor bugs, and I did some improvments
// Mult Rdo module had some serious bugs, but should be ok now.
// BrMultRdo is fixed a lot too.
//
// Revision 1.22  2001/05/31 01:47:23  cholm
// Second rewamp of this directory.  All RDO modules use the common
// BrMultRdoModule, and they both write BrMultRdo Objects.  Also introduced
// ABC for Br[Tile|Si][Parameters|Calibration] since they have a lot in common,
// and the code can be made more efficient this way.
//
// A possible further thing to do, is to make an ABC for the CentModules, and
// the corresponding calibration modules, since they are very VERY similar.
// However, the current module BrMultCentModule is in the way.  Need to talk
// to Steve about that.
//
// Revision 1.21  2001/05/23 12:28:27  cholm
// Uses BrMultModule to obtain vertex
//
// Revision 1.20  2001/05/01 13:00:46  cholm
// Fix "division by zero" problem.
//
// Revision 1.19  2001/04/19 15:14:07  cholm
// Added outlier correction method kHitCorrectopn
//
// Revision 1.18  2001/04/19 00:30:39  sanders
// Removed hard coded BB vertex offset.
//
// Revision 1.17  2001/04/16 02:35:34  sanders
// Changes to have BrSiRdoModule and BrTileRdoModule report the same
// multiplicities as BrRdoModuleMult
//
// Revision 1.16  2001/04/13 03:00:11  sanders
// Further work on reconciling "old" and "new" code.  Added additional checks
// for using TPM1 vertex.
//
// Revision 1.15  2001/04/11 14:01:20  cholm
// Added the method Print for more detailed Printout.
//
// Revision 1.14  2001/04/06 20:00:58  cholm
// Updated to use the per ring conversoj functions as outlined by Steve and
// Hiro in BAN26.
//
// Revision 1.13  2001/04/05 03:05:19  sanders
// Restored "brat complience" to Si and Tile Rdo, BrRdoModuleMult still not complient.
//
// Revision 1.12  2001/04/02 01:18:56  sanders
// BrRdoMult and BrRdoModuleMult reinstated.  These files develop the joint
// Si+Tile multiplicity.  BrRdoModuleMult now gets it parameters from the
// BrTileRdoModule and BrSiRdoModule classes and is a friend of these classes.
//
// Revision 1.11  2001/03/31 16:33:59  bjornhs
// Changed a comment line to remove a compiler problem
//
// Revision 1.10  2001/03/12 19:03:56  cholm
// Revisted the tile classes. Please note, that these are binary incompatible
// with previous classes. These classes now includes (I believe) all that
// Steve and Hiro recently added to thier working code.
//
// Revision 1.7  2001/02/09 18:27:46  cholm
// Removed a left over output from debug phase.
//
// Revision 1.6  2001/02/09 00:11:57  cholm
// Reorganised BrTileRdo and hence also BrTileRdoModule
// Some fixes to BrTileCentModule and BrTileCentCalModule.
//
// Revision 1.5  2001/02/08 15:12:17  cholm
// Fixed a minor bug in BrTileCalibration
// Fleshed out BrtileCent.
// Finalised BrTileCentCalModule.
// Finalised BrTileCentModule.
// Minor bug fixes to BrTileRdoModule.
//
// Revision 1.4  2001/01/31 17:35:17  cholm
// Small corretion in BrTileRdoModule
//
// Revision 1.3  2001/01/31 17:34:24  cholm
// Corrected bug in BrTileParameters when readin from DetectorParameters.txt
//
// Revision 1.2  2001/01/30 23:57:11  cholm
// Corrected some bugs.
//
// Revision 1.1  2001/01/29 20:31:24  cholm
// Split the class BrRdoModuleMult into two different classes, one for the
// silicon, and one for the tiles. Naming corresponds to naming scheme.
//
//














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