Go to the documentation of this file.00001 #include "RecoJets/JetAlgorithms/interface/SubJetAlgorithm.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 #include "fastjet/ClusterSequenceArea.hh"
00005
00006 using namespace std;
00007 using namespace edm;
00008
00009 void SubJetAlgorithm::set_zcut(double z){
00010 zcut_ = z;
00011 }
00012
00013 void SubJetAlgorithm::set_rcut_factor(double r){
00014 rcut_factor_ = r;
00015 }
00016
00017
00018
00019
00020 void SubJetAlgorithm::run( const vector<fastjet::PseudoJet> & cell_particles,
00021 vector<CompoundPseudoJet> & hardjetsOutput ) {
00022
00023
00024
00025 fastjet::FastPrunePlugin PRplugin(*fjJetDefinition_,
00026 *fjJetDefinition_,
00027 zcut_,
00028 rcut_factor_);
00029 fastjet::JetDefinition pjetdef(&PRplugin);
00030
00031
00032
00033
00034
00035 boost::shared_ptr<fastjet::ClusterSequence> fjClusterSeq;
00036 if ( !doAreaFastjet_ ) {
00037 fjClusterSeq =
00038 boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequence( cell_particles,
00039 pjetdef ) );
00040 } else if (voronoiRfact_ <= 0) {
00041 fjClusterSeq =
00042 boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceActiveArea( cell_particles,
00043 pjetdef ,
00044 *fjActiveArea_ ) );
00045 } else {
00046 fjClusterSeq =
00047 boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceVoronoiArea( cell_particles,
00048 pjetdef ,
00049 fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
00050 }
00051
00052 vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(ptMin_);
00053
00054
00055
00056 vector<vector<int> > indices(inclusiveJets.size());
00057
00058 vector<fastjet::PseudoJet>::iterator jetIt = inclusiveJets.begin();
00059 for ( ; jetIt != inclusiveJets.end(); ++jetIt ) {
00060
00061 vector<fastjet::PseudoJet> subjets = fjClusterSeq->exclusive_subjets(*jetIt, nSubjets_);
00062
00063 vector<CompoundPseudoSubJet> subjetsOutput;
00064 std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = subjets.begin(),
00065 itSubJet = itSubJetBegin, itSubJetEnd = subjets.end();
00066 for (; itSubJet != itSubJetEnd; ++itSubJet ){
00067
00068 vector<fastjet::PseudoJet> subjetFastjetConstituents = fjClusterSeq->constituents( *itSubJet );
00069
00070 vector<int> constituents;
00071 vector<fastjet::PseudoJet>::const_iterator fastSubIt = subjetFastjetConstituents.begin(),
00072 transConstBegin = subjetFastjetConstituents.begin(),
00073 transConstEnd = subjetFastjetConstituents.end();
00074 for ( ; fastSubIt != transConstEnd; ++fastSubIt ) {
00075 if (fastSubIt->user_index() >= 0) {
00076 constituents.push_back( fastSubIt->user_index() );
00077 }
00078 }
00079
00080 double subJetArea = (doAreaFastjet_) ?
00081 dynamic_cast<fastjet::ClusterSequenceActiveArea&>(*fjClusterSeq).area(*itSubJet) : 0.0;
00082
00083
00084
00085 subjetsOutput.push_back( CompoundPseudoSubJet( *itSubJet, subJetArea, constituents ) );
00086
00087 }
00088
00089
00090 double fatJetArea = (doAreaFastjet_) ?
00091 dynamic_cast<fastjet::ClusterSequenceActiveArea&>(*fjClusterSeq).area(*jetIt) : 0.0;
00092
00093
00094 hardjetsOutput.push_back( CompoundPseudoJet( *jetIt,fatJetArea,subjetsOutput));
00095 }
00096 }