ATLAS Offline Software
Loading...
Searching...
No Matches
FastJetInterfaceTool.h
Go to the documentation of this file.
1// -*- c++ -*-
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7#ifndef JETREC_FASTJETINTERFACETOOL_H
8#define JETREC_FASTJETINTERFACETOOL_H
9
10//#include "AthenaBaseComps/AthAlgTool.h"
11#include "AsgTools/AsgTool.h"
12
13
14#include "fastjet/ClusterSequence.hh"
15#include "fastjet/JetDefinition.hh"
16#include "fastjet/AreaDefinition.hh"
17
19
20#include <map>
21#include <string>
22#include <vector>
23#include <stdint.h>
24
25#ifndef XAOD_STANDALONE
26// support contribs only in athena for now
27#include "fastjet/SISConePlugin.hh"
28#include "fastjet/CMSIterativeConePlugin.hh"
29#endif
30
31// #include <iostream>
32
34{
36 typedef std::map<std::string,fastjet::JetAlgorithm> algomap_t;
38 typedef std::map<std::string,fastjet::Strategy> strategymap_t;
40 typedef std::map<std::string,fastjet::RecombinationScheme> schememap_t;
42#ifndef XAOD_STANDALONE
43 typedef std::map<std::string,fastjet::SISConePlugin::SplitMergeScale> splitMergeScaleMap_t;
44#else
45 typedef std::map<std::string,int> splitMergeScaleMap_t;
46#endif
48 typedef std::map<std::string,fastjet::AreaType> areamap_t;
49
51 typedef fastjet::PseudoJet fjet_t;
53 typedef std::vector<fjet_t> fjetlist_t;
54
55 // static const algomap_t& getKnownAlgorithms();
56 // static const strategymap_t& getKnownStrategies();
57 // static const schememap_t& getKnownRecombinationSchemes();
58 // static const areamap_t& getKnownAreas();
59
71 template<class D,class M>
72 static bool chkConfig(const std::string& key,D& tag,const M& map)
73 {
74 // std::cout << "Looking for key <" << key << "> in " << map.size() << " entries." << std::endl;
75 typename M::const_iterator fMap(map.find(key));
76 if ( fMap != map.end() ) { tag = (*fMap).second; return true; }
77 else { return false; }
78 }
79
91 template<class D,class M>
92 const std::string& cfgName(D tag,const M& map,const std::string& invalidKey)
93 {
94 typename M::const_iterator fMap(map.begin());
95 typename M::const_iterator lMap(map.end());
96 while ( fMap != lMap && (*fMap).second != tag ) { ++fMap; }
97 return fMap != lMap ? (*fMap).first : invalidKey;
98 }
99}
100
105
106
108 virtual public IFastJetInterfaceTool {
109
110public:
111
113
114
115 FastJetInterfaceTool(const std::string& n);
117 virtual ~FastJetInterfaceTool();
118
123 virtual StatusCode initialize();
124
134 virtual StatusCode execute(const fjetlist_t& inJets,fjetlist_t& outJets);
135
144 virtual fastjet::ClusterSequence* clusterSequence();
145
154 virtual const fastjet::ClusterSequence* clusterSequence() const;
155
167 template<class S>
169 { return dynamic_cast<S*>(this->clusterSequence()); }
170
171
183 template<class S>
184 const S* specificClusterSequence() const
185 { return dynamic_cast<const S*>(this->clusterSequence()); }
186
193 const fastjet::JetDefinition* getJetDefinition() const;
194
201 const fastjet::AreaDefinition* getAreaDefinition() const;
202
209 // cppcheck-suppress returnByReference; do not expose internals
210 const std::string getAreaDefinitionType() const;
211
222
223protected:
224
233
248
250 double m_radius;
258 double m_pTmin;
259
280
288 bool checkConfig(const std::string& key,fastjet::JetAlgorithm& fjalg);
289
297 bool checkConfig(const std::string& key,fastjet::Strategy& fjstr);
298
306 bool checkConfig(const std::string& key,fastjet::RecombinationScheme& fjrs);
307
315#ifndef XAOD_STANDALONE
316 bool checkConfig(const std::string& key,fastjet::SISConePlugin::SplitMergeScale& fjsms);
317#endif
325 bool checkConfig(const std::string& tag,fastjet::AreaType& fjart);
326
333 const std::string& configName(fastjet::JetAlgorithm fjalg);
334
341 const std::string& configName(fastjet::Strategy fjstr);
342
349 const std::string& configName(fastjet::RecombinationScheme fjrs);
350
357#ifndef XAOD_STANDALONE
358 const std::string& configName(fastjet::SISConePlugin::SplitMergeScale fjsms);
359#endif
360
367 const std::string& configName(fastjet::AreaType fjart);
368
369
370private:
371
373 static const std::string m_invalidKeyReference;
374
376 unsigned int m_failedExecCtr;
378 static const unsigned int m_failedExecCtrMax;
379
380
381
382
383protected:
384
386 fastjet::JetAlgorithm m_jetAlgorithm;
388 fastjet::Strategy m_strategy;
390 fastjet::AreaType m_areaType;
392 fastjet::RecombinationScheme m_recombinationScheme;
394#ifndef XAOD_STANDALONE
395 fastjet::SISConePlugin::SplitMergeScale m_SIS_splitMergeScale;
396#endif
398 fastjet::JetDefinition* m_jetDefinition;
400 fastjet::ClusterSequence* m_clusterSequence;
402 fastjet::AreaDefinition* m_areaDefinition;
403
408
409
411 typedef StatusCode
415 typedef StatusCode
417
422
438 StatusCode f_processWithoutArea(const FastJetInterface::fjetlist_t& inJets,
440
455 StatusCode f_processWithArea(const FastJetInterface::fjetlist_t& inJets,
457
464 StatusCode configJetAreas();
498
499
500
502 virtual void updateRandomSeeds();
503
504
505};
506
529
530// inline
531// const FastJetInterface::algomap_t& FastJetInterface::getKnownAlgorithms()
532// { return FastJetInterfaceTool::getKnownAlgorithms(); }
533
534// inline
535// const FastJetInterface::strategymap_t& FastJetInterface::getKnownStrategies()
536// { return FastJetInteraceTool::getKnownStrategies(); }
537
538// inline
539// const FastJetInterface::schememap_t&
540// FastJetInterface::getKnownRecombinationSchemes()
541// { return FastJetInterfaceTool::getKnownRecombinationSchemes(); }
542
543// inline
544// const FastJetInterface::areamap_t& FastJetInterface::getKnownAreas()
545// { return FastJetInterfaceTool::getKnownAreas(); }
546
547#endif
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
double m_exclusiveDcut
Exclusive jet finder property: d cut.
StatusCode(FastJetInterfaceTool::* PROCESSOR)(const FastJetInterface::fjetlist_t &inJets, FastJetInterface::fjetlist_t &outJets)
Processor function pointer type.
std::string m_clusterSequenceType
fastjet property: cluster sequence type
const FastJetInterface::algomap_t & getKnownAlgorithms()
Get known jet algorithms.
std::string m_clusterStrategyType
fastjet property: cluster strategy
fastjet::RecombinationScheme m_recombinationScheme
fastjet tag: recombination scheme
unsigned int m_failedExecCtr
Execution failure counter.
EXTRACTOR m_extractor
Pointer to extractor.
FastJetInterfaceTool(const std::string &n)
Standard AlgTool constructor.
const std::string getAreaDefinitionType() const
Access fastjet area definition type.
static const unsigned int m_failedExecCtrMax
Execution failure counter reporting cut-off.
double m_ghostedMeanKt
fastjet property: average kT of ghosts
uint64_t m_baseRNDSeed
Base seed for random generator. Constructed from m_jetAlgorithm and m_radius.
int m_exclusiveNjets
Exclusive jet finder property: number of jets requested.
std::string m_jetAreaDefinitionType
fastjet property: jet area definition
double m_SIS_protojetPtMin
fastjet property for SIS Cone plugin: protojet min pt
const fastjet::AreaDefinition * getAreaDefinition() const
Access fastjet AreaDefinition.
S * specificClusterSequence()
Access to specific cluster sequence.
virtual StatusCode execute(const fjetlist_t &inJets, fjetlist_t &outJets)
Execute method.
fastjet::AreaDefinition * m_areaDefinition
Pointer to area definition.
double m_ghostedKtScatter
fastjet property: scatter of kT of ghosts
StatusCode configJetAreas()
Configures jet area calculation strategy.
bool m_SIS_doCaching
fastjet property for SIS Cone plugin: do caching
const FastJetInterface::areamap_t & getKnownAreas()
Get known area types.
StatusCode f_extractExclNjets(FastJetInterface::fjetlist_t &outJets)
Exclusive extractor implementation.
static const std::string m_invalidKeyReference
Reference for invalid fastjet tag.
int m_SIS_nPass
fastjet property for SIS Cone plugin: max number of passes
const FastJetInterface::strategymap_t & getKnownStrategies()
Get known jet finding strategies.
const FastJetInterface::splitMergeScaleMap_t & getKnownSplitMergeScales()
Get known split merge scales.
fastjet::Strategy m_strategy
fastjet tag: jet clustering strategy
virtual ~FastJetInterfaceTool()
Base class destructor.
bool checkConfig(const std::string &key, fastjet::JetAlgorithm &fjalg)
Check configuration keyword for jet algorithm.
virtual void updateRandomSeeds()
Set the area rnd seed according to run/event numbers.
const FastJetInterface::schememap_t & getKnownRecombinationSchemes()
Get known recombination schemes.
fastjet::AreaType m_areaType
fastjet tag: jet area type
PROCESSOR m_processor
Pointer to processor.
double m_CMS_seedThreshold
fastjet property for CMS Cone plugin: seed threshold
fastjet::JetDefinition * m_jetDefinition
Pointer to jet definition.
std::string m_jetAlgorithmType
fastjet property: jet algorithm
StatusCode(FastJetInterfaceTool::* EXTRACTOR)(FastJetInterface::fjetlist_t &outJets)
Extractor function pointer type.
bool m_doJetArea
Tool property: flag controls if jet area is calculated.
const fastjet::JetDefinition * getJetDefinition() const
Access fastjet JetDefinition.
bool m_inclusive
Jet algorithm control property: inclusive/exclusive jet finding.
StatusCode f_processWithArea(const FastJetInterface::fjetlist_t &inJets, FastJetInterface::fjetlist_t &outJets)
Processor implementation: cluster sequence with area calculation.
StatusCode f_extractInclusive(FastJetInterface::fjetlist_t &outJets)
Inclusive extractor implementation.
double m_ghostedMinRap
fastjet property: ghost area rapidity limit (min)
double m_ghostedArea
fastjet property: size of area coverded by each ghost
virtual fastjet::ClusterSequence * clusterSequence()
Access fastjet ClusterSequence.
fastjet::ClusterSequence * m_clusterSequence
Pointer to actual cluster sequence.
StatusCode f_processWithoutArea(const FastJetInterface::fjetlist_t &inJets, FastJetInterface::fjetlist_t &outJets)
Processor implementation: cluster sequence without area calculation.
const std::string & configName(fastjet::JetAlgorithm fjalg)
Get keyword for fastjet configuration.
StatusCode f_extractExclDcut(FastJetInterface::fjetlist_t &outJets)
Exclusive extractor implementation.
double m_SIS_splitMergeStopScale
fastjet property for SIS Cone plugin: split merge stop scale
int m_ghostedRepeat
fastjet property: ghost area calculation repeatitions
std::string m_SIS_splitMergeScale_STRING
fastjet property for SIS Cone plugin: split merge scale
const S * specificClusterSequence() const
Access to specific cluster sequence.
virtual StatusCode initialize()
Initialize tool.
fastjet::JetAlgorithm m_jetAlgorithm
fastjet tag: jet algorithm
double m_pTmin
Inclusive jet finder property: pTmin cut.
fastjet::SISConePlugin::SplitMergeScale m_SIS_splitMergeScale
fastjet tag: SIS split merge scale
std::string m_recombinationSchemeType
fastjet property: recombination scheme
double m_ghostedGridScatter
fastjet property: scatter of ghosts on grid
double m_radius
Jet algorithm property: radius/distance parameter.
uint64_t m_userRNDSeed
Seed modifier for random generator : a user property (default to 0, i.e m_baseRNDSeed is unmodified)
double m_ghostedMaxRap
fastjet property: ghost area rapidity limit (max)
double m_SIS_overlapThreshold
fastjet property for SIS Cone plugin: overlap threshold
double m_voronoiEffectiveRfact
fastjet property: Voronoi area effective radius
General tool configuring fastjet functionality.
std::vector< fjet_t > fjetlist_t
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
STL class.
Collection of types for the internal dictionary.
const std::string & cfgName(D tag, const M &map, const std::string &invalidKey)
General configuration keyword finder.
std::map< std::string, fastjet::AreaType > areamap_t
Mapping keyword on fastjet area type tag.
fastjet::PseudoJet fjet_t
fastjet data model for jet
std::map< std::string, fastjet::Strategy > strategymap_t
Mapping keyword on fastjet clustering strategy tag.
std::map< std::string, fastjet::SISConePlugin::SplitMergeScale > splitMergeScaleMap_t
Mapping keyword on fastjet SIS Cone split merge scale tag.
std::map< std::string, fastjet::RecombinationScheme > schememap_t
Mapping keyword on fastjet recombination scheme tag.
static bool chkConfig(const std::string &key, D &tag, const M &map)
General configuration check.
std::map< std::string, fastjet::JetAlgorithm > algomap_t
Mapping keyword on fastjet jet algorithm tag.
std::vector< fjet_t > fjetlist_t
Container for fastjet jets.