CMS 3D CMS Logo

Public Member Functions | Private Attributes

HcalRecHitRecalib Class Reference

#include <CalibCalorimetry/CaloRecalibTools.src/HcalRecHitRecalib.cc>

Inheritance diagram for HcalRecHitRecalib:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 HcalRecHitRecalib (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &) override
 ~HcalRecHitRecalib ()

Private Attributes

edm::InputTag hbheLabel_
std::string hcalfile_
std::string hcalfileinpath_
edm::InputTag hfLabel_
edm::InputTag hoLabel_
CaloMiscalibMapHcal mapHcal_
std::string RecalibHBHEHits_
std::string RecalibHFHits_
std::string RecalibHOHits_
double refactor_
double refactor_mean_

Detailed Description

Description: Producer to miscalibrate (calibrated) Hcal RecHit

Implementation: <Notes on="" implementation>="">

Definition at line 38 of file HcalRecHitRecalib.h.


Constructor & Destructor Documentation

HcalRecHitRecalib::HcalRecHitRecalib ( const edm::ParameterSet iConfig) [explicit]

Definition at line 12 of file HcalRecHitRecalib.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hbheLabel_, hcalfile_, hcalfileinpath_, hfLabel_, hoLabel_, RecalibHBHEHits_, RecalibHFHits_, RecalibHOHits_, refactor_, refactor_mean_, and AlCaHLTBitMon_QueryRunRegistry::string.

{

  hbheLabel_ = iConfig.getParameter<edm::InputTag>("hbheInput");
  hoLabel_ = iConfig.getParameter<edm::InputTag>("hoInput");
  hfLabel_ = iConfig.getParameter<edm::InputTag>("hfInput");

//   HBHEHitsProducer_ = iConfig.getParameter< std::string > ("HBHERecHitsProducer");
//   HOHitsProducer_ = iConfig.getParameter< std::string > ("HERecHitsProducer");
//   HFHitsProducer_ = iConfig.getParameter< std::string > ("HERecHitsProducer");
//   HBHEHits_ = iConfig.getParameter< std::string > ("HBHEHitCollection");
//   HFHits_ = iConfig.getParameter< std::string > ("HFHitCollection");
//   HOHits_ = iConfig.getParameter< std::string > ("HOHitCollection");

  RecalibHBHEHits_ = iConfig.getParameter< std::string > ("RecalibHBHEHitCollection");
  RecalibHFHits_ = iConfig.getParameter< std::string > ("RecalibHFHitCollection");
  RecalibHOHits_ = iConfig.getParameter< std::string > ("RecalibHOHitCollection");

  refactor_ = iConfig.getUntrackedParameter<double> ("Refactor",(double)1);
  refactor_mean_ = iConfig.getUntrackedParameter<double> ("Refactor_mean",(double)1);

  //register your products
  produces< HBHERecHitCollection >(RecalibHBHEHits_);
  produces< HFRecHitCollection >(RecalibHFHits_);
  produces< HORecHitCollection >(RecalibHOHits_);

  // here read them from xml (particular to HCAL)

  hcalfileinpath_=iConfig.getUntrackedParameter<std::string> ("fileNameHcal","");
  edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);

  hcalfile_=hcalfiletmp.fullPath();
}
HcalRecHitRecalib::~HcalRecHitRecalib ( )

Definition at line 47 of file HcalRecHitRecalib.cc.

{
 

}

Member Function Documentation

void HcalRecHitRecalib::beginRun ( const edm::Run ,
const edm::EventSetup iSetup 
) [override, virtual]
void HcalRecHitRecalib::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [override, virtual]

Implements edm::EDProducer.

Definition at line 68 of file HcalRecHitRecalib.cc.

References CaloMiscalibMapHcal::get(), edm::Event::getByLabel(), hbheLabel_, egHLT::errCodes::HBHERecHits, hfLabel_, egHLT::errCodes::HFRecHits, hoLabel_, LogDebug, mapHcal_, edm::Event::put(), RecalibHBHEHits_, RecalibHFHits_, RecalibHOHits_, refactor_, and refactor_mean_.

{
  using namespace edm;
  using namespace std;

  Handle<HBHERecHitCollection> HBHERecHitsHandle;
  Handle<HFRecHitCollection> HFRecHitsHandle;
  Handle<HORecHitCollection> HORecHitsHandle;

  const HBHERecHitCollection*  HBHERecHits = 0;
  const HFRecHitCollection*  HFRecHits = 0;
  const HORecHitCollection*  HORecHits = 0;

  iEvent.getByLabel(hbheLabel_,HBHERecHitsHandle);
  if (!HBHERecHitsHandle.isValid()) {
    LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
  } else {
    HBHERecHits = HBHERecHitsHandle.product(); // get a ptr to the product
  }

  iEvent.getByLabel(hoLabel_,HORecHitsHandle);
  if (!HORecHitsHandle.isValid()) {
    LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
  } else {
    HORecHits = HORecHitsHandle.product(); // get a ptr to the product
  }

  iEvent.getByLabel(hfLabel_,HFRecHitsHandle);
  if (!HFRecHitsHandle.isValid()) {
    LogDebug("") << "HcalREcHitRecalib: Error! can't get product!" << std::endl;
  } else {
    HFRecHits = HFRecHitsHandle.product(); // get a ptr to the product
  }


//     iEvent.getByLabel(HBHEHitsProducer_,HBHEHits_,HBHERecHitsHandle);
//     HBHERecHits = HBHERecHitsHandle.product(); // get a ptr to the product

//     iEvent.getByLabel(HFHitsProducer_,HFHits_,HFRecHitsHandle);
//     HFRecHits = HFRecHitsHandle.product(); // get a ptr to the product

//     iEvent.getByLabel(HOHitsProducer_,HOHits_,HORecHitsHandle);
//     HORecHits = HORecHitsHandle.product(); // get a ptr to the product


  //Create empty output collections
  std::auto_ptr< HBHERecHitCollection > RecalibHBHERecHitCollection( new HBHERecHitCollection );
  std::auto_ptr< HFRecHitCollection > RecalibHFRecHitCollection( new HFRecHitCollection );
  std::auto_ptr< HORecHitCollection > RecalibHORecHitCollection( new HORecHitCollection );

  // Intercalib constants
  //  edm::ESHandle<EcalIntercalibConstants> pIcal;
  //  iSetup.get<EcalIntercalibConstantsRcd>().get(pIcal);
  //  const EcalIntercalibConstants* ical = pIcal.product();

  if(HBHERecHits)
    {

       //loop on all EcalRecHits (barrel)
      HBHERecHitCollection::const_iterator itHBHE;
      for (itHBHE=HBHERecHits->begin(); itHBHE!=HBHERecHits->end(); itHBHE++) {
        
        // find intercalib constant for this cell

        //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
        //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;

        //      if( icalit!=ical->getMap().end() ){
        //        icalconst = icalit->second;
          // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;

        //      } else {
        //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
        //              ;
        //          }
          
        float icalconst=(mapHcal_.get().find(itHBHE->id().rawId()))->second;
          // make the rechit with rescaled energy and put in the output collection

        icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
        HBHERecHit aHit(itHBHE->id(),itHBHE->energy()*icalconst,itHBHE->time());
        
        RecalibHBHERecHitCollection->push_back( aHit);
      }
    }

  if(HFRecHits)
    {

       //loop on all EcalRecHits (barrel)
      HFRecHitCollection::const_iterator itHF;
      for (itHF=HFRecHits->begin(); itHF!=HFRecHits->end(); itHF++) {
        
        // find intercalib constant for this cell

        //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
        //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;

        //      if( icalit!=ical->getMap().end() ){
        //        icalconst = icalit->second;
          // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;

        //      } else {
        //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
        //              ;
        //          }
          
          // make the rechit with rescaled energy and put in the output collection
        
        float icalconst=(mapHcal_.get().find(itHF->id().rawId()))->second;
        icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
        HFRecHit aHit(itHF->id(),itHF->energy()*icalconst,itHF->time());
        
        RecalibHFRecHitCollection->push_back( aHit);
      }
    }

  if(HORecHits)
    {

       //loop on all EcalRecHits (barrel)
      HORecHitCollection::const_iterator itHO;
      for (itHO=HORecHits->begin(); itHO!=HORecHits->end(); itHO++) {
        
        // find intercalib constant for this cell

        //      EcalIntercalibConstants::EcalIntercalibConstantMap::const_iterator icalit=ical->getMap().find(itb->id().rawId());
        //      EcalIntercalibConstants::EcalIntercalibConstant icalconst;

        //      if( icalit!=ical->getMap().end() ){
        //        icalconst = icalit->second;
          // edm::LogDebug("EcalRecHitMiscalib") << "Found intercalib for xtal " << EBDetId(itb->id()) << " " << icalconst ;

        //      } else {
        //        edm::LogError("EcalRecHitMiscalib") << "No intercalib const found for xtal " << EBDetId(itb->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
        //              ;
        //          }
          
          // make the rechit with rescaled energy and put in the output collection

        float icalconst=(mapHcal_.get().find(itHO->id().rawId()))->second;
        icalconst=refactor_mean_+(icalconst-refactor_mean_)*refactor_; //apply additional scaling factor (works if gaussian)    
        HORecHit aHit(itHO->id(),itHO->energy()*icalconst,itHO->time());
          
          RecalibHORecHitCollection->push_back( aHit);
      }
    }


  //Put Recalibrated rechit in the event
  iEvent.put( RecalibHBHERecHitCollection, RecalibHBHEHits_);
  iEvent.put( RecalibHFRecHitCollection, RecalibHFHits_);
  iEvent.put( RecalibHORecHitCollection, RecalibHOHits_);
}

Member Data Documentation

Definition at line 48 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

std::string HcalRecHitRecalib::hcalfile_ [private]

Definition at line 55 of file HcalRecHitRecalib.h.

Referenced by beginRun(), and HcalRecHitRecalib().

std::string HcalRecHitRecalib::hcalfileinpath_ [private]

Definition at line 56 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib().

Definition at line 50 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

Definition at line 49 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

Definition at line 58 of file HcalRecHitRecalib.h.

Referenced by beginRun(), and produce().

std::string HcalRecHitRecalib::RecalibHBHEHits_ [private]

Definition at line 51 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

std::string HcalRecHitRecalib::RecalibHFHits_ [private]

Definition at line 52 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

std::string HcalRecHitRecalib::RecalibHOHits_ [private]

Definition at line 53 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

double HcalRecHitRecalib::refactor_ [private]

Definition at line 59 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().

Definition at line 60 of file HcalRecHitRecalib.h.

Referenced by HcalRecHitRecalib(), and produce().