CMS 3D CMS Logo

CambridgeAlgorithmWrapper Class Reference

CambridgeAlgorithmWrapper is the Wrapper subclass which runs the CambridgeAlgorithm of FastJet for jetfinding. More...

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

Inheritance diagram for CambridgeAlgorithmWrapper:

FastJetBaseWrapper

List of all members.

Public Member Functions

 CambridgeAlgorithmWrapper (const edm::ParameterSet &fConfig)
virtual ~CambridgeAlgorithmWrapper ()


Detailed Description

CambridgeAlgorithmWrapper is the Wrapper subclass which runs the CambridgeAlgorithm of FastJet for jetfinding.

The FastJet package, written by Matteo Cacciari and Gavin Salam, provides a fast implementation of the longitudinally invariant kt and longitudinally invariant inclusive Cambridge/Aachen jet finders. More information can be found at: https://parthe.lpthe.jussieu.fr/~salam/fastjet/

Authors:
Andreas Oehler, University Karlsruhe (TH) and Dorian Kcira, Institut de Physique Nucleaire Departement de Physique Universite Catholique de Louvain have written the CambridgeAlgorithmWrapper class which uses the above mentioned package within the Framework of CMSSW
Version:
1st Version Nov. 6 2006

2nd redesigned by F.Ratnikov (UMd) Aug. 1, 2007

Definition at line 33 of file CambridgeAlgorithmWrapper.h.


Constructor & Destructor Documentation

CambridgeAlgorithmWrapper::CambridgeAlgorithmWrapper ( const edm::ParameterSet fConfig  ) 

Definition at line 18 of file CambridgeAlgorithmWrapper.cc.

References lat::endl(), edm::ParameterSet::getParameter(), and FastJetBaseWrapper::mJetDefinition.

00019   : FastJetBaseWrapper (fConfig)
00020 {
00021   //configuring algorithm 
00022   double rParam = fConfig.getParameter<double> ("FJ_ktRParam");
00023   //choosing search-strategy:
00024   fastjet::Strategy fjStrategy = fastjet::plugin_strategy;
00025   std::string strategy = fConfig.getParameter<std::string>("Strategy");
00026   if (strategy == "Best") fjStrategy = fastjet::Best; // Chooses best Strategy for every event, depending on N and ktRParam
00027   else if (strategy == "N2Plain") fjStrategy = fastjet::N2Plain;  // N2Plain is best for N<50
00028   else if (strategy == "N2Tiled") fjStrategy = fastjet::N2Tiled;  // N2Tiled is best for 50<N<400
00029   else if (strategy == "N2MinHeapTiled") fjStrategy = fastjet::N2MinHeapTiled;  // N2MinHeapTiles is best for 400<N<15000
00030   else if (strategy == "NlnN") fjStrategy = fastjet::NlnN;  // NlnN is best for N>15000
00031   else if (strategy == "NlnNCam") fjStrategy = fastjet::NlnNCam;  // NlnNCam is best for N>6000
00032   else  edm::LogError("FastJetDefinition") << "CambridgeAlgorithmWrapper-> Unknown strategy: " 
00033                                            << strategy
00034                                            << ". Known strategies: Best N2Plain N2Tiled N2MinHeapTiled NlnN" << std::endl;
00035   //additional strategies are possible, but not documented in the manual as they are experimental,
00036   //they are also not used by the "Best" method. Note: "NlnNCam" only works with 
00037   //the cambridge_algorithm and does not link against CGAL, "NlnN" is only 
00038   //available if fastjet is linked against CGAL.
00039   //The above given numbers of N are for ktRaram=1.0, for other ktRParam the numbers would be 
00040   //different.
00041   
00042   mJetDefinition = new fastjet::JetDefinition (fastjet::cambridge_algorithm, rParam, fjStrategy);
00043 }

CambridgeAlgorithmWrapper::~CambridgeAlgorithmWrapper (  )  [virtual]

Definition at line 45 of file CambridgeAlgorithmWrapper.cc.

00045 {}


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