BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//__________________________________________________________________
//                                                                  
//	BrChkvRdoModule                                              
//                                                                   
//      Module class that takes digitized data (BrChkvDig)          
//      and calibration parameters to reconstruct Chkv hits         
//                                                                   
//__________________________________________________________________

//
// $Id: BrChkvRdoModule.cxx,v 1.11 2002/02/14 15:13:24 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/02/14 15:13:24 $
// $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//

#ifndef WIN32
#include <iostream>
#include <iomanip>
#else
#include <iostream.h>
#include <iomanip.h>
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TNtuple
#include "TNtuple.h"
#endif
#ifndef BRAT_BrChkvRdoModule
#include "BrChkvRdoModule.h"
#endif
#ifndef BRAT_BrGeometryDbManager
#include "BrGeometryDbManager.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrC1Dig
#include "BrC1Dig.h"
#endif
#ifndef BRAT_BrRichDig
#include "BrRichDig.h"
#endif
#ifndef BRAT_BrChkvRdo
#include "BrChkvRdo.h"
#endif

//___________________________________________________________
ClassImp(BrChkvRdoModule);

//___________________________________________________________
 BrChkvRdoModule::BrChkvRdoModule()
  : BrModule()
{
  // Default Constructor. (DO NOT USE)
  fChkvParams   = 0;
  fCalParams   = 0;
  fChkvVol      = 0;

  SetThreshold(0); 
  SetNoPedWidthCut(0);
}

//___________________________________________________________
 BrChkvRdoModule::BrChkvRdoModule(Text_t *Name,Char_t *Title)
  : BrModule(Name,Title)
{
  fChkvParams   = 0;
  fCalParams   = 0;
  fChkvVol      = 0;

  if (strcmp(GetName(),"C1")==0) {
    SetThreshold(0); 
    SetNoPedWidthCut(0);
  }
  if (strcmp(GetName(),"RICH")==0) {
    SetThreshold(0); 
    SetNoPedWidthCut(3.0);
  }
}

//___________________________________________________________
 BrChkvRdoModule::~BrChkvRdoModule()
{
  // destructor
}

//___________________________________________________________
 void BrChkvRdoModule::DefineHistograms()
{
  // define histograms here
  TDirectory* saveDir = gDirectory;
  
  TDirectory* dir = gDirectory->mkdir(Form("%sRdo",GetName()));
  dir->cd();
  
  // histos
  Int_t noTubes = fChkvParams->GetNoTubes();

  fNHits           = new TH1F("nHits","Number of calibrated hits", 
			      noTubes, 0.5, noTubes+0.5);
  fAdcVsTube       = new TH2F("adcVsTube","Raw Adc versus Tube",
			      noTubes, 0.5, noTubes+0.5, 2000, 0, 4000);
  fCalAdcVsTube    = new TH2F("calAdcVsTube","Calibrated Adc versus Tube",
			      noTubes, 0.5, noTubes+0.5, 220, -1, 10);
  fCalAdcVsTubeAcc = new TH2F("calAdcVsTubeAcc","Accepted Calibrated Adc versus Tube",
			      noTubes, 0.5, noTubes+0.5, 550, -1, 10);  

  fNHits->SetXTitle("Number of accepted hits per event");
  fAdcVsTube->SetXTitle("Tube number");
  fAdcVsTube->SetYTitle("Raw Adc");
  fCalAdcVsTube->SetXTitle("Tube number");
  fCalAdcVsTube->SetYTitle("Calibrated Adc");
  fCalAdcVsTubeAcc->SetXTitle("Tube number");
  fCalAdcVsTubeAcc->SetYTitle("Calibrated Adc");

  // list of histograms
  dir->mkdir("Tubes");
  dir->cd("Tubes");      
  
  fTubeAdc = new TH1F*[noTubes];
  for (Int_t t=1; t<=noTubes; t++) {
    fTubeAdc[t-1] = new TH1F(Form("adcTube%02d",t),
			     Form("Raw Adc Tube %02d",t),
			     2000,0,4000);
    fTubeAdc[t-1]->SetXTitle("Raw Adc");
  }

  fTubeCalAdc = new TH1F*[noTubes];
  for (Int_t t=1; t<=noTubes; t++) {
    fTubeCalAdc[t-1] = new TH1F(Form("calAdcTube%02d",t),
				Form("Calibrated Adc Tube %02d",t),
				420,-1,20);
    fTubeCalAdc[t-1]->SetXTitle("Calibrated Adc");
  }

  gDirectory = saveDir;  
}

//___________________________________________________________
 void BrChkvRdoModule::Init(){
  //-------------------
  // initialize module
  //-------------------

  if(DebugLevel() > 1)
    cout << "Entering Init of BrChkvRdoModule" << endl;
  
  // get important information from managers

  BrCalibrationManager *parMan =
    BrCalibrationManager::Instance();
  
  if (!parMan) {
    Failure("Init", 
	    "Couldn't instantiate calibration manager");
    return;
  }  
  
  // get the calibration element (created in main)
  fCalParams =
    (BrChkvCalibration*)parMan->Register("BrChkvCalibration", GetName());
  
  if (!fCalParams) {
    Failure("Init", 
	    "Couldn't get calibration parameters for %s", GetName());
    return;
  }
  
  fCalParams->Use("pedestal");
  fCalParams->Use("pedestalWidth");
  fCalParams->Use("adcGain");  
  
  BrGeometryDbManager  *geoDb = BrGeometryDbManager::Instance();
  BrParameterDbManager *parDb = BrParameterDbManager::Instance();
  
  if (!geoDb) {
    Failure("Init ", "Couldn't instantiate geometry manager");
    return;
  }
  
  if (!parDb) {
    Failure("Init ", "Couldn't instantiate  detector parameter manager");
    return;
  }
  
  fChkvParams = 
    (BrChkvParameters*)parDb->GetDetectorParameters("BrChkvParameters",
						    (Char_t*)GetName());
  
  fChkvVol  = (BrDetectorVolume*)geoDb->
    GetDetectorVolume("BrDetectorVolume", GetName());
  if(!fChkvVol)
    Abort("Init","No Cherenkov Detector Volume : %s", GetName());

}

//____________________________________________________________________
 void BrChkvRdoModule::Begin()
{
  SetState(kBegin);

  // check if calibration revision exist
  if (!fCalParams->RevisionExists("pedestal")) {
    Fatal("Begin","Could not find pedestal calibration revision!");
    return;
  }
  if (!fCalParams->RevisionExists("pedestal")) {
    Fatal("Begin","Could not find adcgain calibration revision!");
    return;
  }
}

//____________________________________________________________________
 void BrChkvRdoModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // -----------------------------------------------------------------
  // Event method to be called once per event.
  // 1- get raw cherenkov data 
  // 2- get calibration, set adc and position
  // 3- store recontructed hits in rdo object in out node
  // -----------------------------------------------------------------

  // The event method is splittet up in two (identical) sections
  // one for C1 and one for thr RICH - it could not be done in
  // one section since BrC1Dig and BrRichDig does not have a 
  // common baseclass (which is a bit annoying!).

  SetState(kEvent);
  
  BrC1Dig* c1Dig = 0;
  BrRichDig* richDig = 0;


  // C1 section ++++++++++++++++++++++++++++++++++++++++++++++++++++
  if (strcmp(GetName(),"C1")==0) {
    c1Dig = (BrC1Dig*) inNode->GetObject("DigC1");

    if (!c1Dig) {
      Warning("Event","DigC1 not found!");
      return;
    }  

    // prepare output rdo object
    BrChkvRdo* c1Rdo = new BrChkvRdo(Form("ChkvRdo %s", GetName()),
				     Form("ChkvRdo %s", GetName()));
    
    // loop over C1 tubes
    Int_t nHits = 0;
    for(Int_t t = 1; t <= fChkvParams->GetNoTubes(); t++){
      
      //get calibration parameters
      Double_t ped = fCalParams->GetPedestal(t);
      Double_t pedW = fCalParams->GetPedestalWidth(t);
      Double_t gain = fCalParams->GetAdcGain(t);

      if ((ped==BrChkvCalibration::kCalException)||
	  (pedW==BrChkvCalibration::kCalException)||
	  (gain==BrChkvCalibration::kCalException)) continue;

      //calibrate!
      Double_t adc    = c1Dig->GetAdc(t);
      Double_t calAdc = (adc - ped)/gain;
      
      if(fHistOn) {
	fTubeAdc[t-1]->Fill(adc-ped);
	fAdcVsTube->Fill(t,adc);
      	fCalAdcVsTube->Fill(t,calAdc);
      }

      // Accept hit??? 
      if ( calAdc < fThreshold ) continue;
      if ( (adc - ped) < (fNoPedWidthCut * pedW) ) continue;
      nHits++;

      if(fHistOn) {
	fTubeCalAdc[t-1]->Fill(calAdc);
	fCalAdcVsTubeAcc->Fill(t,calAdc);
      }
      
      c1Rdo->AddHit(t, calAdc);
    } // end of c1 tubes loop

    if(fHistOn)
      fNHits->Fill(nHits);

    if (c1Rdo->GetEntries())
      outNode->AddObject(c1Rdo);
    else
      delete c1Rdo;
  } 
  // end of C1 section +++++++++++++++++++++++++++++++++++++++++++++++    
  
  // RICH section ++++++++++++++++++++++++++++++++++++++++++++++++++++
  if (strcmp(GetName(),"RICH")==0) {
    richDig = (BrRichDig*) inNode->GetObject("DigRICH");
    
    if (!richDig) {
      Warning("Event","DigRICH not found!");
      return;
    }

    // prepare output rdo object
    BrChkvRdo* richRdo = new BrChkvRdo(Form("ChkvRdo %s", GetName()),
				       Form("ChkvRdo %s", GetName()));
    
    // loop over RICH tubes
    Int_t nHits = 0;
    for(Int_t t = 1; t <= fChkvParams->GetNoTubes(); t++){
      
      
      //get calibration parameters
      Double_t ped = fCalParams->GetPedestal(t);
      Double_t pedW = fCalParams->GetPedestalWidth(t);
      Double_t gain = fCalParams->GetAdcGain(t);

      if ((ped==BrChkvCalibration::kCalException)||
	  (pedW==BrChkvCalibration::kCalException)||
	  (gain==BrChkvCalibration::kCalException)) continue;

      //calibrate!
      Double_t adc    = richDig->GetAdc(t);
      Double_t calAdc = (adc - ped) / gain;
      
      if(fHistOn) {
	fTubeAdc[t-1]->Fill(adc-ped);
	fAdcVsTube   ->Fill(t,adc);
      	fCalAdcVsTube->Fill(t,calAdc);
      }

      // Accept hit??? 
      //if ( calAdc < fThreshold ) continue;
      // only cut in n times pedestalwidth over pedestal
      if ( (adc - ped) < (fNoPedWidthCut * pedW) ) continue; 
      nHits++;
      
      if(fHistOn) {
	fTubeCalAdc[t-1]->Fill(calAdc);
	fCalAdcVsTubeAcc->Fill(t,calAdc);
      }
      
      richRdo->AddHit(t, calAdc);      
    } // end of rich tubes loop
    
    if(fHistOn)
      fNHits->Fill(nHits);
    
    if (richRdo->GetEntries())
      outNode->AddObject(richRdo);
    else
      delete richRdo;

  } 
  // end of RICH section +++++++++++++++++++++++++++++++++++++++++++++
}

//__________________________________________________________________
 void BrChkvRdoModule::Print(Option_t* option) const
{
  // Standard information printout. 
  // 
  // Options: See BrModule::Print
  // 
  TString opt(option);
  opt.ToLower();
  
  BrModule::Print(option);
  if (opt.Contains("d"))
    cout << endl
	 << "  Chkv Reconstruction Module" << endl
         << "  Original author: Claus E. Jorgensen" << endl
         << "  Revisted by:     $Author: videbaek $" << endl
         << "  Revision date:   $Date: 2002/02/14 15:13:24 $"   << endl
         << "  Revision number: $Revision: 1.11 $ " << endl
         << endl
         << "*************************************************" << endl;
}


//
// $Log: BrChkvRdoModule.cxx,v $
// Revision 1.11  2002/02/14 15:13:24  videbaek
// Added protection in Init for no volume
//
// Revision 1.10  2002/02/13 11:50:53  ekman
// Introduced adc gain calibration of the RICH
//
// Revision 1.9  2001/11/21 16:49:07  ekman
// Changed adc histograms (now with pedestal subtraction).
//
// Revision 1.8  2001/10/17 16:34:30  ekman
// Introduced histograms for each tube and deleted the ntuple with tube
// info.
//
// Revision 1.7  2001/10/17 00:22:32  ekman
// Introduced a new cut method to select hits - now one can select by cutting
// relative to cal adc = 1 (threshold) or requiring adc > ped + n* pedW. Updated the histograms.
//
// Revision 1.6  2001/10/08 11:29:08  cholm
// Changed to use new DB access classes
//
// Revision 1.5  2001/09/24 13:11:29  ekman
// Requires (2.5 * pedestal) width for valid hit in RICH tubes.
//
// Revision 1.4  2001/09/04 10:47:16  ekman
// Init: Only use adcGain if the detector is C1. Gain calibrations for
// the Rich are not available yet!
//
// Revision 1.3  2001/08/27 16:43:50  ekman
// Changed histogram ranges (again).
//
// Revision 1.2  2001/08/21 19:43:31  ekman
// Added adc gain calibration. Changed histogram ranges.
//
// Revision 1.1  2001/08/09 16:01:33  ekman
// Added the classes BrChkvRdo and BrChkvRcoModule.
//
//

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