  | 
  
    ATLAS Offline Software
    
   | 
 
 
 
 
Go to the documentation of this file.
   35                                           const std::string& 
name,
 
   36                                           const IInterface* 
parent ) :
 
   39   declareInterface< DerivationFramework::IAugmentationTool >(
this);
 
   73                  << 
" Assuming it's an auxdata member");
 
   78   return StatusCode::SUCCESS;
 
   88   return StatusCode::SUCCESS;
 
   94 #define COPY_AND_SORT_CONTAINER( CONTAINERTYPE )                                       \ 
   95 else if ( evtStore()->contains<CONTAINERTYPE>( m_inCollKey.value() ) ) {               \ 
   96   ATH_MSG_DEBUG("Trying to copy, sort, and record container of type "#CONTAINERTYPE ); \ 
   97   const CONTAINERTYPE* inCont;                                                         \ 
   98   ATH_CHECK( evtStore()->retrieve( inCont, m_inCollKey.value() ) );                    \ 
   99   CONTAINERTYPE* outCont = new CONTAINERTYPE( SG::VIEW_ELEMENTS );                     \ 
  100   *outCont = *inCont;                                                                  \ 
  101   ATH_CHECK( evtStore()->record ( outCont, m_outCollKey.value() ) );                   \ 
  102   ATH_CHECK( this->doSort(outCont) );                                                  \ 
  107 #define OVERWRITE_AND_SORT_CONTAINER( CONTAINERTYPE )                                                \ 
  108 else if ( evtStore()->contains<CONTAINERTYPE>( m_inCollKey.value() ) ) {                             \ 
  109   ATH_MSG_DEBUG("Trying to copy, sort, and overwrite container of type "#CONTAINERTYPE );            \ 
  110   const CONTAINERTYPE* inCont;                                                                       \ 
  111   ATH_CHECK( evtStore()->retrieve( inCont, m_inCollKey.value() ) );                                  \ 
  112   ConstDataVector<CONTAINERTYPE>* outCont = new ConstDataVector<CONTAINERTYPE>( SG::VIEW_ELEMENTS ); \ 
  113   for ( const CONTAINERTYPE::base_value_type* inPart : *inCont ){                                    \ 
  114     outCont->push_back(inPart);                                                                      \ 
  116   ATH_CHECK( evtStore()->overwrite( outCont, m_inCollKey.value() ) );                                \ 
  117   ATH_CHECK( this->doSortConst<CONTAINERTYPE>(outCont) );                                            \ 
  133                   << 
"Trying to retrieve a non-const version of the 'InputContainer'...");
 
  137       ATH_MSG_DEBUG(
"We couldn't retrieve a non-const version of the input container... try const.");
 
  156         ATH_MSG_ERROR(
"Couln't find the provided intput container in store gate for later overwriting");
 
  157         return StatusCode::FAILURE;
 
  162     ATH_MSG_DEBUG(
"Got a non-empty 'OutputCollection' property. " 
  163                   << 
"Trying to retrieve a const version of the 'InputContainer'...");
 
  181       ATH_MSG_ERROR(
"Couln't find the provided intput container in store gate");
 
  182       return StatusCode::FAILURE;
 
  187   return StatusCode::SUCCESS;
 
  196     return StatusCode::FAILURE;
 
  205   else if ( std::abs(
m_sortID) == 2 ) {
 
  210   else if ( std::abs(
m_sortID) == 3 ) {
 
  215   else if ( std::abs(
m_sortID) == 4 ) {
 
  220   else if ( std::abs(
m_sortID) == 5 ) {
 
  225   else if ( std::abs(
m_sortID) == 6 ) {
 
  230   else if ( std::abs(
m_sortID) == 7 ) {
 
  236   return StatusCode::SUCCESS;
 
  243   const double a = partA->
pt();
 
  244   const double b = partB->
pt();
 
  252   const double a = partA->
eta();
 
  253   const double b = partB->
eta();
 
  261   const double a = partA->
phi();
 
  262   const double b = partB->
phi();
 
  270   const double a = partA->
m();
 
  271   const double b = partB->
m();
 
  279   const double a = partA->
e();
 
  280   const double b = partB->
e();
 
  
def retrieve(aClass, aKey=None)
 
virtual ~ParticleSortingTool()
Destructor:
 
bool compareMass(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's mass.
 
bool compareAuxData(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare an auxdata member of the particle.
 
DataVector adapter that acts like it holds const pointers.
 
StringProperty m_outCollKey
The name of the output container (with SG::VIEW_ELEMENTS) with the sorted copy of input objects.
 
int m_sortID
Internal identifier for the type of sorting.
 
StringProperty m_sortVar
Define by what parameter to sort (default: 'pt')
 
virtual double rapidity() const =0
The true rapidity (y) of the particle.
 
#define COPY_AND_SORT_CONTAINER(CONTAINERTYPE)
 
bool comparePt(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's pt.
 
std::atomic< unsigned long > m_nEventsProcessed
Internal event counter.
 
Class providing the definition of the 4-vector interface.
 
bool compareRapidity(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's rapidity.
 
virtual StatusCode initialize() override
Athena algtool's initialize.
 
bool comparePhi(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's phi.
 
ServiceHandle< StoreGateSvc > & evtStore()
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
 
bool compareEta(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's eta.
 
ParticleSortingTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
 
::StatusCode StatusCode
StatusCode definition for legacy code.
 
bool compareEnergy(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's energy.
 
BooleanProperty m_sortDescending
Define if the container should be sorted in a descending order (default=true)
 
virtual double pt() const =0
The transverse momentum ( ) of the particle.
 
StringProperty m_inCollKey
Input container name.
 
bool compareDouble(double a, double b) const
Method to compare two doubles.
 
void sort()
Sort the container.
 
virtual double eta() const =0
The pseudorapidity ( ) of the particle.
 
StatusCode doSort(xAOD::IParticleContainer *cont) const
Helper method that implements the call to the right sort function.
 
virtual StatusCode finalize() override
Athena algtool's finalize.
 
virtual double phi() const =0
The azimuthal angle ( ) of the particle.
 
virtual StatusCode addBranches(const EventContext &ctx) const final override
Implement the method from the ISkimmingTool interface.
 
virtual double e() const =0
The total energy of the particle.
 
virtual double m() const =0
The invariant mass of the particle.
 
#define OVERWRITE_AND_SORT_CONTAINER(CONTAINERTYPE)