CMS 3D CMS Logo

FastJetBaseWrapper Class Reference

#include <RecoJets/JetAlgorithms/interface/FastJetBaseWrapper.h>

Inheritance diagram for FastJetBaseWrapper:

CambridgeAlgorithmWrapper CDFMidpointAlgorithmWrapper KtJetAlgorithmWrapper SISConeAlgorithmWrapper

List of all members.

Public Member Functions

 FastJetBaseWrapper (const edm::ParameterSet &fConfig)
void run (const JetReco::InputCollection &fInput, JetReco::OutputCollection *fOutput)
virtual ~FastJetBaseWrapper ()

Protected Attributes

fastjet::JetDefinition * mJetDefinition

Private Attributes

fastjet::GhostedAreaSpec * mActiveArea
double mJetPtMin


Detailed Description

Definition at line 13 of file FastJetBaseWrapper.h.


Constructor & Destructor Documentation

FastJetBaseWrapper::FastJetBaseWrapper ( const edm::ParameterSet fConfig  ) 

Definition at line 24 of file FastJetBaseWrapper.cc.

References edm::ParameterSet::getParameter(), mActiveArea, and mJetPtMin.

00025   : mJetDefinition (0), 
00026     mActiveArea (0)
00027 {
00028   mJetPtMin = fConfig.getParameter<double> ("jetPtMin");
00029   if (fConfig.getParameter<std::string> ("UE_Subtraction") == "yes") {            // accept pilup subtraction parameters
00030     double ghostEtaMax = fConfig.getParameter<double> ("Ghost_EtaMax");   //default Ghost_EtaMax should be 6
00031     int activeAreaRepeats = fConfig.getParameter<int> ("Active_Area_Repeats");   //default Active_Area_Repeats 5
00032     double ghostArea = fConfig.getParameter<double> ("GhostArea");   //default GhostArea 0.01
00033     mActiveArea = new fastjet::ActiveAreaSpec (ghostEtaMax, activeAreaRepeats, ghostArea);
00034   }
00035 }

FastJetBaseWrapper::~FastJetBaseWrapper (  )  [virtual]

Definition at line 38 of file FastJetBaseWrapper.cc.

References mActiveArea, and mJetDefinition.

00038                                        {
00039   delete mJetDefinition;
00040   delete mActiveArea;
00041 }


Member Function Documentation

void FastJetBaseWrapper::run ( const JetReco::InputCollection fInput,
JetReco::OutputCollection fOutput 
)

Definition at line 43 of file FastJetBaseWrapper.cc.

References c, i, iggi_31X_cfg::input, metsig::jet, pfTauBenchmarkGeneric_cfi::jets, mActiveArea, mJetDefinition, mJetPtMin, and p4.

Referenced by cms::KtPilupSubtractionJetProducer::runAlgorithm().

00043                                                                                                    {
00044   if (!fOutput) return;
00045   fOutput->clear();
00046   
00047   // convert inputs
00048   std::vector<fastjet::PseudoJet> fjInputs;
00049   fjInputs.reserve (fInput.size());
00050   
00051   JetReco::InputCollection::const_iterator input = fInput.begin();
00052   for (unsigned i = 0; i < fInput.size(); ++i) {
00053     const JetReco::InputItem& c = fInput[i];
00054     fjInputs.push_back (fastjet::PseudoJet (c->px(),c->py(),c->pz(),c->energy()));
00055     fjInputs.back().set_user_index(i);
00056   }
00057    
00058    // create an object that represents your choice of jet finder and 
00059    // the associated parameters
00060    // run the jet clustering with the above jet definition
00061 
00062    // here we need to keep both pointers, as "area" interfaces are missing in base class
00063    fastjet::ClusterSequenceActiveArea* clusterSequenceWithArea = 0;
00064    fastjet::ClusterSequenceWithArea* clusterSequence = 0;
00065    if (mActiveArea) {
00066      clusterSequenceWithArea = new fastjet::ClusterSequenceActiveArea (fjInputs, *mJetDefinition, *mActiveArea);
00067      clusterSequence = clusterSequenceWithArea;
00068    }
00069    else {
00070      clusterSequence = new fastjet::ClusterSequenceWithArea (fjInputs, *mJetDefinition);
00071    }
00072    // retrieve jets for selected mode
00073    std::vector<fastjet::PseudoJet> jets = clusterSequence->inclusive_jets (mJetPtMin);
00074 
00075    // get PU pt
00076    double median_Pt_Per_Area = clusterSequenceWithArea ? clusterSequenceWithArea->pt_per_unit_area() : 0.;
00077 
00078    // process found jets
00079    for (std::vector<fastjet::PseudoJet>::const_iterator jet=jets.begin(); jet!=jets.end();++jet) {
00080      // jet itself
00081      double px=jet->px();
00082      double py=jet->py();
00083      double pz=jet->pz();
00084      double E=jet->E();
00085      double jetArea=clusterSequence->area(*jet);
00086      double pu=0.;
00087      // PU subtraction
00088      if (clusterSequenceWithArea) {
00089        fastjet::PseudoJet pu_p4 = median_Pt_Per_Area * clusterSequenceWithArea->area_4vector(*jet);
00090        pu = pu_p4.E();
00091        if (pu_p4.perp2() >= jet->perp2() || pu_p4.E() >= jet->E()) { // if the correction is too large, set the jet to zero
00092          px = py = pz = E = 0.;
00093        } 
00094        else {   // otherwise do an E-scheme subtraction
00095          px -= pu_p4.px();
00096          py -= pu_p4.py();
00097          pz -= pu_p4.pz();
00098          E -= pu_p4.E();
00099        }
00100      }
00101      math::XYZTLorentzVector p4(px,py,pz,E);
00102      // constituents
00103      std::vector<fastjet::PseudoJet> fastjet_constituents = clusterSequence->constituents(*jet);
00104      JetReco::InputCollection jetConstituents; 
00105      jetConstituents.reserve (fastjet_constituents.size());
00106      for (std::vector<fastjet::PseudoJet>::const_iterator itConst=fastjet_constituents.begin();
00107           itConst!=fastjet_constituents.end();itConst++){
00108        jetConstituents.push_back(fInput[(*itConst).user_index()]);
00109      }
00110      // Build ProtoJet
00111      fOutput->push_back(ProtoJet(p4,jetConstituents));
00112      fOutput->back().setJetArea (jetArea);
00113      fOutput->back().setPileup (pu);
00114    }
00115    // cleanup
00116    if (clusterSequenceWithArea) delete clusterSequenceWithArea;
00117    else delete clusterSequence; // sigh... No plymorphism in fastjet
00118 }


Member Data Documentation

fastjet::GhostedAreaSpec* FastJetBaseWrapper::mActiveArea [private]

Definition at line 21 of file FastJetBaseWrapper.h.

Referenced by FastJetBaseWrapper(), run(), and ~FastJetBaseWrapper().

fastjet::JetDefinition* FastJetBaseWrapper::mJetDefinition [protected]

Definition at line 19 of file FastJetBaseWrapper.h.

Referenced by CambridgeAlgorithmWrapper::CambridgeAlgorithmWrapper(), CDFMidpointAlgorithmWrapper::CDFMidpointAlgorithmWrapper(), KtJetAlgorithmWrapper::KtJetAlgorithmWrapper(), run(), SISConeAlgorithmWrapper::SISConeAlgorithmWrapper(), and ~FastJetBaseWrapper().

double FastJetBaseWrapper::mJetPtMin [private]

Definition at line 22 of file FastJetBaseWrapper.h.

Referenced by FastJetBaseWrapper(), and run().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:20:42 2009 for CMSSW by  doxygen 1.5.4