CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Attributes

HLTRHemisphere Class Reference

#include <HLTRHemisphere.h>

Inheritance diagram for HLTRHemisphere:
HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 HLTRHemisphere (const edm::ParameterSet &)
 ~HLTRHemisphere ()

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Private Attributes

bool accNJJets_
edm::InputTag inputTag_
double max_Eta_
int max_NJ_
double min_Jet_Pt_
bool saveTag_

Detailed Description

Definition at line 16 of file HLTRHemisphere.h.


Constructor & Destructor Documentation

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

Definition at line 30 of file HLTRHemisphere.cc.

References accNJJets_, edm::InputTag::encode(), inputTag_, LogDebug, max_Eta_, max_NJ_, and min_Jet_Pt_.

                                                             :
  inputTag_    (iConfig.getParameter<edm::InputTag>("inputTag")),
  min_Jet_Pt_  (iConfig.getParameter<double>       ("minJetPt" )),
  max_Eta_     (iConfig.getParameter<double>       ("maxEta" )),
  max_NJ_      (iConfig.getParameter<int>          ("maxNJ" )),
  accNJJets_   (iConfig.getParameter<bool>         ("acceptNJ" ))
{
   LogDebug("") << "Input/minJetPt/maxEta/maxNJ/acceptNJ : "
                << inputTag_.encode() << " "
                << min_Jet_Pt_ << "/"
                << max_Eta_ << "/"
                << max_NJ_ << "/"
                << accNJJets_ << ".";

   //register your products
   produces<std::vector<math::XYZTLorentzVector> >();
}
HLTRHemisphere::~HLTRHemisphere ( )

Definition at line 48 of file HLTRHemisphere.cc.

{
}

Member Function Documentation

void HLTRHemisphere::fillDescriptions ( edm::ConfigurationDescriptions descriptions) [static]

Reimplemented from edm::EDFilter.

Definition at line 53 of file HLTRHemisphere.cc.

References edm::ParameterSetDescription::add(), and edm::ConfigurationDescriptions::add().

                                                                           {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
  desc.add<double>("minJetPt",30.0);
  desc.add<double>("maxEta",3.0);
  desc.add<int>("maxNJ",7);
  desc.add<bool>("acceptNJ",true);
  descriptions.add("hltRHemisphere",desc);
}
bool HLTRHemisphere::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements HLTFilter.

Definition at line 69 of file HLTRHemisphere.cc.

References accNJJets_, prof2calltree::count, edm::Event::getByLabel(), i, inputTag_, analyzePatCleaning_cfg::jets, max_Eta_, max_NJ_, min_Jet_Pt_, n, edm::Event::put(), and dt_offlineAnalysis_common_cff::reco.

{
   using namespace std;
   using namespace edm;
   using namespace reco;
   using namespace math;

   // get hold of collection of objects
   Handle<CaloJetCollection> jets;
   iEvent.getByLabel (inputTag_,jets);

   // The output Collection
   std::auto_ptr<vector<math::XYZTLorentzVector> > Hemispheres(new vector<math::XYZTLorentzVector> );

   // look at all objects, check cuts and add to filter object
   int n(0);
   reco::CaloJetCollection JETS;
   CaloJetCollection::const_iterator i ( jets->begin() );
   for (unsigned int i=0; i<jets->size(); i++) {
     if(fabs(jets->at(i).eta()) < max_Eta_ && jets->at(i).pt() >= min_Jet_Pt_){
       JETS.push_back(jets->at(i));
       n++;
     }
   }

  if(n<2){
    return false; //need at least 2 jets to build the hemispheres
  }

  if(n>max_NJ_ && max_NJ_!=-1){
    iEvent.put(Hemispheres);
    return accNJJets_; // 
  }
   int N_comb(1); // compute the number of combinations of jets possible
  for(unsigned int i = 0; i < JETS.size(); i++){
    N_comb *= 2;                
  }
  //Make the hemispheres
  XYZTLorentzVector j1,j2;
  double M_min = 9999999999.0;
  double dHT_min = 99999999.0;
  int j_count;
  for(int i=0;i<N_comb;i++){       
    XYZTLorentzVector j_temp1, j_temp2;
    int itemp = i;
    j_count = N_comb/2;
    int count = 0;
    while(j_count > 0){
      if(itemp/j_count == 1){
        j_temp1 += JETS.at(count).p4();
      } else {
        j_temp2 += JETS.at(count).p4();
      }
      itemp -= j_count*(itemp/j_count);
      j_count /= 2;
      count++;
    }
    double M_temp = j_temp1.M2()+j_temp2.M2();
    if(M_temp < M_min){
      M_min = M_temp;
      j1= j_temp1;
      j2= j_temp2;
    }
    double dHT_temp = fabs(j_temp1.E()-j_temp2.E());
    if(dHT_temp < dHT_min){
      dHT_min = dHT_temp;
    }
  }

  Hemispheres->push_back(j1);
  Hemispheres->push_back(j2);

  iEvent.put(Hemispheres);
  return true;
}

Member Data Documentation

Definition at line 31 of file HLTRHemisphere.h.

Referenced by filter(), and HLTRHemisphere().

Definition at line 26 of file HLTRHemisphere.h.

Referenced by filter(), and HLTRHemisphere().

double HLTRHemisphere::max_Eta_ [private]

Definition at line 29 of file HLTRHemisphere.h.

Referenced by filter(), and HLTRHemisphere().

int HLTRHemisphere::max_NJ_ [private]

Definition at line 30 of file HLTRHemisphere.h.

Referenced by filter(), and HLTRHemisphere().

double HLTRHemisphere::min_Jet_Pt_ [private]

Definition at line 28 of file HLTRHemisphere.h.

Referenced by filter(), and HLTRHemisphere().

bool HLTRHemisphere::saveTag_ [private]

Definition at line 27 of file HLTRHemisphere.h.