CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Attributes

HLTRFilter Class Reference

#include <HLTRFilter.h>

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

List of all members.

Public Member Functions

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

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Private Attributes

bool accept_NJ_
bool DoRPrime_
edm::InputTag inputMetTag_
edm::InputTag inputTag_
double min_MR_
double min_R_
bool saveTag_

Detailed Description

Definition at line 16 of file HLTRFilter.h.


Constructor & Destructor Documentation

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

Definition at line 29 of file HLTRFilter.cc.

References accept_NJ_, DoRPrime_, edm::InputTag::encode(), inputMetTag_, inputTag_, LogDebug, min_MR_, and min_R_.

                                                     :
  inputTag_    (iConfig.getParameter<edm::InputTag>("inputTag")),
  inputMetTag_ (iConfig.getParameter<edm::InputTag>("inputMetTag")),
  min_R_       (iConfig.getParameter<double>       ("minR"   )),
  min_MR_      (iConfig.getParameter<double>       ("minMR"   )),
  DoRPrime_    (iConfig.getParameter<bool>       ("doRPrime"   )),
  accept_NJ_    (iConfig.getParameter<bool>       ("acceptNJ"   ))

{
   LogDebug("") << "Inputs/minR/minMR/doRPrime/acceptNJ : "
                << inputTag_.encode() << " "
                << inputMetTag_.encode() << " "
                << min_R_ << " "
                << min_MR_ << " "
                << DoRPrime_ << " "
                << accept_NJ_ << ".";

   //register your products
   produces<trigger::TriggerFilterObjectWithRefs>();
}
HLTRFilter::~HLTRFilter ( )

Definition at line 50 of file HLTRFilter.cc.

{
}

Member Function Documentation

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

Reimplemented from edm::EDFilter.

Definition at line 55 of file HLTRFilter.cc.

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

                                                                       {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltRHemisphere"));
  desc.add<edm::InputTag>("inputMetTag",edm::InputTag("hltMet"));
  desc.add<double>("minR",0.3);
  desc.add<double>("minMR",100.0);
  desc.add<bool>("doRPrime",false);
  desc.add<bool>("acceptNJ",true);
  descriptions.add("hltRFilter",desc);
}
bool HLTRFilter::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements HLTFilter.

Definition at line 72 of file HLTRFilter.cc.

References accept_NJ_, edm::Event::getByLabel(), inputMetTag_, inputTag_, CaloMET_cfi::met, min_MR_, min_R_, module(), path(), phi, edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, mathSSE::sqrt(), and cond::rpcobtemp::temp.

{
   using namespace std;
   using namespace edm;
   using namespace reco;
   using namespace trigger;

   // The filter object
   auto_ptr<TriggerFilterObjectWithRefs>
     filterobject (new TriggerFilterObjectWithRefs(path(),module()));

   // get hold of collection of objects
   Handle< vector<math::XYZTLorentzVector> > hemispheres;
   iEvent.getByLabel (inputTag_,hemispheres);

   // get hold of the MET Collection
   Handle<CaloMETCollection> inputMet;
   iEvent.getByLabel(inputMetTag_,inputMet);  

   iEvent.put(filterobject);

   // check the the input collections are available
   if (not hemispheres.isValid() or not inputMet.isValid())
     return false;

   if(hemispheres->size() ==0){  // the Hemisphere Maker will produce an empty collection of hemispheres if the number of jets in the
     return accept_NJ_;   // event is greater than the maximum number of jets
   }

   //***********************************
   //Calculate R or R prime

   TLorentzVector j1R(hemispheres->at(0).x(),hemispheres->at(0).y(),hemispheres->at(0).z(),hemispheres->at(0).t());
   TLorentzVector j2R(hemispheres->at(1).x(),hemispheres->at(1).y(),hemispheres->at(1).z(),hemispheres->at(1).t());
   TLorentzVector j1Rp(hemispheres->at(2).x(),hemispheres->at(2).y(),hemispheres->at(2).z(),hemispheres->at(2).t());
   TLorentzVector j2Rp(hemispheres->at(3).x(),hemispheres->at(3).y(),hemispheres->at(3).z(),hemispheres->at(3).t());

   if(j1R.Pt() > 0.1) { // some good R combination was found

     j1R.SetPtEtaPhiM(j1R.Pt(),j1R.Eta(),j1R.Phi(),0.0);
     j2R.SetPtEtaPhiM(j2R.Pt(),j2R.Eta(),j2R.Phi(),0.0);
  
     if(j2R.Pt() > j1R.Pt()){
       TLorentzVector temp = j1R;
       j1R = j2R;
       j2R = temp;
     }
  
     double num = j1R.P()-j2R.P();
     double den = j1R.Pz()-j2R.Pz();
     if(fabs(num)<fabs(den)) {

       //now we can calculate MTR
       TVector3 met;
       met.SetPtEtaPhi((inputMet->front()).pt(),0.0,(inputMet->front()).phi());
       double MTR = sqrt(0.5*(met.Mag()*(j1R.Pt()+j2R.Pt()) - met.Dot(j1R.Vect()+j2R.Vect())));

       //calculate MR
       double MR=0;
       double temp = (j1R.P()*j2R.Pz()-j2R.P()*j1R.Pz())*(j1R.P()*j2R.Pz()-j2R.P()*j1R.Pz());
       temp /= (j1R.Pz()-j2R.Pz())*(j1R.Pz()-j2R.Pz())-(j1R.P()-j2R.P())*(j1R.P()-j2R.P());    
       MR = 2.*sqrt(temp);
       if(MR>=min_MR_ && float(MTR)/float(MR)>=min_R_) return true;
     }
   }

   if(j1Rp.Pt() > 0.1) { // some good R' combination was found  
     
     j1Rp.SetPtEtaPhiM(j1Rp.Pt(),j1Rp.Eta(),j1Rp.Phi(),0.0);
     j2Rp.SetPtEtaPhiM(j2Rp.Pt(),j2Rp.Eta(),j2Rp.Phi(),0.0);

     if(j2Rp.Pt() > j1Rp.Pt()){
       TLorentzVector temp = j1Rp;
       j1Rp = j2Rp;
       j2Rp = temp;
     }

     double num = j1Rp.P()-j2Rp.P();
     double den = j1Rp.Pz()-j2Rp.Pz();
     if(fabs(num)>fabs(den)) {
       //now we can calculate MTR

       TVector3 met;
       met.SetPtEtaPhi((inputMet->front()).pt(),0.0,(inputMet->front()).phi());
       double MTR = sqrt(0.5*(met.Mag()*(j1Rp.Pt()+j2Rp.Pt()) - met.Dot(j1Rp.Vect()+j2Rp.Vect())));

       double jaP = j1Rp.Pt()*j1Rp.Pt() +j1Rp.Pz()*j2Rp.Pz()-j1Rp.P()*j2Rp.P();
       double jbP = j2Rp.Pt()*j2Rp.Pt() +j1Rp.Pz()*j2Rp.Pz()-j1Rp.P()*j2Rp.P();
       jbP *= -1.;
       double den = sqrt((j1Rp.P()-j2Rp.P())*(j1Rp.P()-j2Rp.P())-(j1Rp.Pz()-j2Rp.Pz())*(j1Rp.Pz()-j2Rp.Pz()));
       
       jaP /= den;
       jbP /= den;
    
       double temp = jaP*met.Dot(j2Rp.Vect())/met.Mag() + jbP*met.Dot(j1Rp.Vect())/met.Mag();
       temp = temp*temp;
    
       den = (met.Dot(j1Rp.Vect()+j2Rp.Vect())/met.Mag())*(met.Dot(j1Rp.Vect()+j2Rp.Vect())/met.Mag())-(jaP-jbP)*(jaP-jbP);
    
       if(den <= 0.0) return false;
       
       temp /= den;
       temp = 2.*sqrt(temp);
       
       double bR = (jaP-jbP)/(met.Dot(j1Rp.Vect()+j2Rp.Vect())/met.Mag());
       double gR = 1./sqrt(1.-bR*bR);
    
       temp *= gR;

       double MRp = temp;

       if(MRp>=min_MR_ && float(MTR)/float(MRp)>=min_R_) return true;
     }
   }

   // filter decision
   return false;
}

Member Data Documentation

bool HLTRFilter::accept_NJ_ [private]

Definition at line 32 of file HLTRFilter.h.

Referenced by filter(), and HLTRFilter().

bool HLTRFilter::DoRPrime_ [private]

Definition at line 31 of file HLTRFilter.h.

Referenced by HLTRFilter().

Definition at line 27 of file HLTRFilter.h.

Referenced by filter(), and HLTRFilter().

Definition at line 26 of file HLTRFilter.h.

Referenced by filter(), and HLTRFilter().

double HLTRFilter::min_MR_ [private]

Definition at line 30 of file HLTRFilter.h.

Referenced by filter(), and HLTRFilter().

double HLTRFilter::min_R_ [private]

Definition at line 29 of file HLTRFilter.h.

Referenced by filter(), and HLTRFilter().

bool HLTRFilter::saveTag_ [private]

Definition at line 28 of file HLTRFilter.h.