CMS 3D CMS Logo

Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes

HiL1Subtractor Class Reference

#include <RecoHI/HiJetAlgos/plugins/HiL1Subtractor.cc>

Inheritance diagram for HiL1Subtractor:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 HiL1Subtractor (const edm::ParameterSet &)
 ~HiL1Subtractor ()

Protected Attributes

std::string jetType_
std::string rhoTag_

Private Member Functions

virtual void beginJob ()
virtual void endJob ()
virtual void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

edm::InputTag src_

Detailed Description

Description:

Implementation:

Definition at line 42 of file HiL1Subtractor.h.


Constructor & Destructor Documentation

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

Definition at line 23 of file HiL1Subtractor.cc.

References Exception, and jetType_.

                                                             :
  src_       ( iConfig.getParameter<edm::InputTag>("src") ),
  jetType_       (iConfig.getParameter<std::string>("jetType") ),
  rhoTag_       (iConfig.getParameter<std::string>("rhoTag") )

{

  if(jetType_ == "CaloJet"){
    produces<reco::CaloJetCollection >();
  }
  else if(jetType_ == "PFJet"){
    produces<reco::PFJetCollection >();
  }
  else if(jetType_ == "GenJet"){
     produces<reco::GenJetCollection >();
  }
  else{
    throw cms::Exception("InvalidInput") << "invalid jet type in HiL1Subtractor\n";
  }

}
HiL1Subtractor::~HiL1Subtractor ( )

Definition at line 46 of file HiL1Subtractor.cc.

{
}

Member Function Documentation

void HiL1Subtractor::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 254 of file HiL1Subtractor.cc.

{
}
void HiL1Subtractor::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 260 of file HiL1Subtractor.cc.

                       {
}
void HiL1Subtractor::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 57 of file HiL1Subtractor.cc.

References reco::LeafCandidate::et(), reco::LeafCandidate::eta(), edm::Event::getByLabel(), j, metsig::jet, reco::Jet::jetArea(), analyzePatCleaning_cfg::jets, jetType_, medianPtkt, edm::Event::put(), rho, rhoTag_, reco::Jet::scaleEnergy(), reco::Jet::setPileup(), and src_.

{

  // get the input jet collection and create output jet collection

  // right now, identical loop for calo and PF jets, should template
   if(jetType_ == "GenJet"){
      std::auto_ptr<reco::GenJetCollection> jets( new reco::GenJetCollection);
      edm::Handle< edm::View<reco::GenJet> > h_jets;
      iEvent.getByLabel( src_, h_jets );

      // Grab appropriate rho, hard coded for the moment                                                                                      
      edm::Handle<std::vector<double> > rs;
      iEvent.getByLabel(edm::InputTag(rhoTag_,"rhos"),rs);
      //iEvent.getByLabel(edm::InputTag("ic5CaloJets","rhos"),rs);                                                                            
      int rsize = rs->size();

      for(int j = 0; j < rsize; j++){
         double medianpt=rs->at(j);
         medianPtkt[j]=medianpt;
      }

      // loop over the jets                                                                                                                   
      int jetsize = h_jets->size();
      for(int ijet = 0; ijet < jetsize; ++ijet){

         reco::GenJet jet = ((*h_jets)[ijet]);

         double jet_eta = jet.eta();
         double jet_et = jet.et();

         //std::cout<<" pre-subtracted jet_et "<<jet_et<<std::endl;                                                                            

         if(fabs(jet_eta)<=3){

            double rho=-999;

            if (jet_eta<-2.5 && jet_eta>-3.5)rho=medianPtkt[2];
            if (jet_eta<-1.5 && jet_eta>-2.5)rho=medianPtkt[3];
            if (jet_eta<-0.5 && jet_eta>-1.5)rho=medianPtkt[4];
            if (jet_eta<0.5 && jet_eta>-0.5)rho=medianPtkt[5];
            if (jet_eta<1.5 && jet_eta>0.5)rho=medianPtkt[6];
            if (jet_eta<2.5 && jet_eta>1.5)rho=medianPtkt[7];
            if (jet_eta<3.5 && jet_eta>2.5)rho=medianPtkt[8];

            double jet_area = jet.jetArea();

            double CorrFactor =0.;
            if(rho*jet_area<jet_et) CorrFactor = 1.0 - rho*jet_area/jet_et;
            jet.scaleEnergy( CorrFactor );
            jet.setPileup(rho*jet_area);
            
            //std::cout<<"  correction factor "<<1.0 - rho*jet_area/jet_et<<std::endl;                                                          
         }

         //std::cout<<" subtracted jet_et "<<jet.et()<<std::endl;                                                                              
         jets->push_back(jet);


      }
      iEvent.put(jets);

   }else if(jetType_ == "CaloJet"){
    std::auto_ptr<reco::CaloJetCollection> jets( new reco::CaloJetCollection);
    edm::Handle< edm::View<reco::CaloJet> > h_jets;
    iEvent.getByLabel( src_, h_jets );

    // Grab appropriate rho, hard coded for the moment
    edm::Handle<std::vector<double> > rs;
    iEvent.getByLabel(edm::InputTag(rhoTag_,"rhos"),rs);
    //iEvent.getByLabel(edm::InputTag("ic5CaloJets","rhos"),rs);


    int rsize = rs->size();

    for(int j = 0; j < rsize; j++){
      double medianpt=rs->at(j);
      medianPtkt[j]=medianpt;
    }

    // loop over the jets

    int jetsize = h_jets->size();

    for(int ijet = 0; ijet < jetsize; ++ijet){

      reco::CaloJet jet = ((*h_jets)[ijet]);

      double jet_eta = jet.eta();
      double jet_et = jet.et();

      //std::cout<<" pre-subtracted jet_et "<<jet_et<<std::endl;

      if(fabs(jet_eta)<=3){

        double rho=-999;

        if (jet_eta<-2.5 && jet_eta>-3.5)rho=medianPtkt[2];
        if (jet_eta<-1.5 && jet_eta>-2.5)rho=medianPtkt[3];
        if (jet_eta<-0.5 && jet_eta>-1.5)rho=medianPtkt[4];
        if (jet_eta<0.5 && jet_eta>-0.5)rho=medianPtkt[5];
        if (jet_eta<1.5 && jet_eta>0.5)rho=medianPtkt[6];
        if (jet_eta<2.5 && jet_eta>1.5)rho=medianPtkt[7];
        if (jet_eta<3.5 && jet_eta>2.5)rho=medianPtkt[8];

        double jet_area = jet.jetArea();

        double CorrFactor =0.;
        if(rho*jet_area<jet_et) CorrFactor = 1.0 - rho*jet_area/jet_et;
        jet.scaleEnergy( CorrFactor );
        jet.setPileup(rho*jet_area);

        //std::cout<<"  correction factor "<<1.0 - rho*jet_area/jet_et<<std::endl;
      }

      //std::cout<<" subtracted jet_et "<<jet.et()<<std::endl;

      jets->push_back(jet);


    }
    iEvent.put(jets);

  }
  else if(jetType_ == "PFJet"){
    std::auto_ptr<reco::PFJetCollection> jets( new reco::PFJetCollection);
    edm::Handle< edm::View<reco::PFJet> > h_jets;
    iEvent.getByLabel( src_, h_jets );

    // Grab appropriate rho, hard coded for the moment
    edm::Handle<std::vector<double> > rs;
    iEvent.getByLabel(edm::InputTag(rhoTag_,"rhos"),rs);
    //iEvent.getByLabel(edm::InputTag("ic5CaloJets","rhos"),rs);


    int rsize = rs->size();

    for(int j = 0; j < rsize; j++){
      double medianpt=rs->at(j);
      medianPtkt[j]=medianpt;
    }

    // loop over the jets

    int jetsize = h_jets->size();

    for(int ijet = 0; ijet < jetsize; ++ijet){

      reco::PFJet jet = ((*h_jets)[ijet]);

      double jet_eta = jet.eta();
      double jet_et = jet.et();

      //std::cout<<" pre-subtracted jet_et "<<jet_et<<std::endl;

      if(fabs(jet_eta)<=3){

        double rho=-999;

        if (jet_eta<-2.5 && jet_eta>-3.5)rho=medianPtkt[2];
        if (jet_eta<-1.5 && jet_eta>-2.5)rho=medianPtkt[3];
        if (jet_eta<-0.5 && jet_eta>-1.5)rho=medianPtkt[4];
        if (jet_eta<0.5 && jet_eta>-0.5)rho=medianPtkt[5];
        if (jet_eta<1.5 && jet_eta>0.5)rho=medianPtkt[6];
        if (jet_eta<2.5 && jet_eta>1.5)rho=medianPtkt[7];
        if (jet_eta<3.5 && jet_eta>2.5)rho=medianPtkt[8];

        double jet_area = jet.jetArea();

        double CorrFactor =0.;
        if(rho*jet_area<jet_et) CorrFactor = 1.0 - rho*jet_area/jet_et;
        jet.scaleEnergy( CorrFactor );
        jet.setPileup(rho*jet_area);

        //std::cout<<"  correction factor "<<1.0 - rho*jet_area/jet_et<<std::endl;
      }

      //std::cout<<" subtracted jet_et "<<jet.et()<<std::endl;

      jets->push_back(jet);


    }
    iEvent.put(jets);

  }







}

Member Data Documentation

std::string HiL1Subtractor::jetType_ [protected]

Definition at line 65 of file HiL1Subtractor.h.

Referenced by HiL1Subtractor(), and produce().

std::string HiL1Subtractor::rhoTag_ [protected]

Definition at line 66 of file HiL1Subtractor.h.

Referenced by produce().

Definition at line 62 of file HiL1Subtractor.h.

Referenced by produce().