FLAVOR_T

FlavorHistory

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115
#ifndef HepMCCandidate_FlavorHistory_h
#define HepMCCandidate_FlavorHistory_h

/** \class reco::FlavorHistory
 *
 * Stores information about the flavor history of a parton
 *
 * \author: Stephen Mrenna (FNAL), Salvatore Rappoccio (JHU)
 *
 */

// -------------------------------------------------------------
// Identify the ancestry of the Quark
//
//
// Matrix Element:
//    Status 3 parent with precisely 2 "grandparents" that
//    is outside of the "initial" section (0-5) that has the
//    same ID as the status 2 parton in question.
//
// Flavor excitation:
//    If we find only one outgoing parton.
//
// Gluon splitting:
//    Parent is a quark of a different flavor than the parton
//    in question, or a gluon.
//    Can come from either ISR or FSR.
//
// True decay:
//    Decays from a resonance like top, Higgs, etc.
// -------------------------------------------------------------

#include "DataFormats/Common/interface/Ptr.h"
#include "DataFormats/Common/interface/OwnVector.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/Candidate/interface/Candidate.h"
#include "DataFormats/Candidate/interface/CandidateFwd.h"
#include "DataFormats/Candidate/interface/ShallowClonePtrCandidate.h"

#include <fstream>

namespace reco {

  class FlavorHistory {
  public:
    enum FLAVOR_T {
      FLAVOR_NULL = 0,  // No flavor, unset
      FLAVOR_GS,        // gluon splitting
      FLAVOR_EXC,       // flavor excitation
      FLAVOR_ME,        // matrix element
      FLAVOR_DECAY,     // flavor decay
      N_FLAVOR_TYPES
    };  // total number

    static const int gluonId = 21;
    static const int tQuarkId = 6;
    static const int bQuarkId = 5;
    static const int cQuarkId = 4;

    FlavorHistory();
    FlavorHistory(FLAVOR_T flavorSource,
                  reco::CandidatePtr const& parton,
                  reco::CandidatePtr const& progenitor,
                  reco::CandidatePtr const& sister,
                  reco::ShallowClonePtrCandidate const& matchedJet,
                  reco::ShallowClonePtrCandidate const& sisterJet);
    FlavorHistory(FLAVOR_T flavorSource,
                  edm::Handle<edm::View<reco::Candidate> > h_partons,
                  int iparton,
                  int iprogenitor,
                  int isister,
                  reco::ShallowClonePtrCandidate const& matchedJet,
                  reco::ShallowClonePtrCandidate const& sisterJet);
    FlavorHistory(FLAVOR_T flavorSource,
                  edm::Handle<reco::CandidateCollection> h_partons,
                  int iparton,
                  int iprogenitor,
                  int isister,
                  reco::ShallowClonePtrCandidate const& matchedJet,
                  reco::ShallowClonePtrCandidate const& sisterJet);
    ~FlavorHistory() {}

    // Accessors
    FLAVOR_T flavorSource() const { return flavorSource_; }
    bool hasParton() const { return parton_.isNonnull(); }
    bool hasSister() const { return sister_.isNonnull(); }
    bool hasProgenitor() const { return progenitor_.isNonnull(); }
    bool hasMatchedJet() const { return matchedJet_.masterClonePtr().isNonnull(); }
    bool hasSisterJet() const { return sisterJet_.masterClonePtr().isNonnull(); }
    const reco::CandidatePtr& parton() const { return parton_; }
    const reco::CandidatePtr& sister() const { return sister_; }
    const reco::CandidatePtr& progenitor() const { return progenitor_; }
    const reco::ShallowClonePtrCandidate& matchedJet() const { return matchedJet_; }
    const reco::ShallowClonePtrCandidate& sisterJet() const { return sisterJet_; }

    // Operators for sorting and keys
    bool operator<(FlavorHistory const& right) const { return parton_.key() < right.parton_.key(); }
    bool operator>(FlavorHistory const& right) const { return parton_.key() > right.parton_.key(); }
    bool operator==(FlavorHistory const& right) const { return parton_.key() == right.parton_.key(); }

  protected:
    FLAVOR_T flavorSource_;
    reco::CandidatePtr parton_;
    reco::CandidatePtr progenitor_;
    reco::CandidatePtr sister_;
    reco::ShallowClonePtrCandidate matchedJet_;
    reco::ShallowClonePtrCandidate sisterJet_;
  };

}  // namespace reco

std::ostream& operator<<(std::ostream& out, reco::FlavorHistory const& cand);

#endif