ParticleState

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 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219
#ifndef Candidate_ParticleState_h
#define Candidate_ParticleState_h
/** \class reco::ParticleState
 *
 * Base class describing a generic reconstructed particle
 * its main subclass is Candidate
 *
 * \author Luca Lista, INFN
 *
 *
 */

#include <Rtypes.h>

#include "DataFormats/Math/interface/Point3D.h"
#include "DataFormats/Math/interface/Vector3D.h"
#include "DataFormats/Math/interface/PtEtaPhiMass.h"
#include "DataFormats/GeometryVector/interface/GlobalVector.h"
#include "DataFormats/Math/interface/LorentzVector.h"

namespace reco {

  class ParticleState {
  public:
    /// electric charge type
    typedef int Charge;
    /// Lorentz vector
    typedef math::XYZTLorentzVector LorentzVector;
    /// Lorentz vector
    typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
    /// point in the space
    typedef math::XYZPoint Point;
    /// point in the space
    typedef math::XYZVector Vector;
    /// default constructor
    ParticleState() : vertex_(0, 0, 0), qx3_(0), pdgId_(0), status_(0) {}

    /// constructor from values
    ParticleState(Charge q,
                  const PtEtaPhiMass& p4,
                  const Point& vertex = Point(0, 0, 0),
                  int pdgId = 0,
                  int status = 0,
                  bool integerCharge = true)
        : vertex_(vertex),
          p4Polar_(p4.pt(), p4.eta(), p4.phi(), p4.mass()),
          p4Cartesian_(p4Polar_),
          qx3_(integerCharge ? q * 3 : q),
          pdgId_(pdgId),
          status_(status) {}

    /// constructor from values
    ParticleState(Charge q,
                  const LorentzVector& p4,
                  const Point& vertex = Point(0, 0, 0),
                  int pdgId = 0,
                  int status = 0,
                  bool integerCharge = true)
        : vertex_(vertex),
          p4Polar_(p4),
          p4Cartesian_(p4),
          qx3_(integerCharge ? q * 3 : q),
          pdgId_(pdgId),
          status_(status) {}

    /// constructor from values
    ParticleState(Charge q,
                  const PolarLorentzVector& p4,
                  const Point& vertex = Point(0, 0, 0),
                  int pdgId = 0,
                  int status = 0,
                  bool integerCharge = true)
        : vertex_(vertex),
          p4Polar_(p4),
          p4Cartesian_(p4),
          qx3_(integerCharge ? q * 3 : q),
          pdgId_(pdgId),
          status_(status) {}

    ParticleState(Charge q,
                  const GlobalVector& p3,
                  float iEnergy,
                  float imass,
                  const Point& vertex = Point(0, 0, 0),
                  int pdgId = 0,
                  int status = 0,
                  bool integerCharge = true)
        : vertex_(vertex),
          p4Polar_(p3.perp(), p3.eta(), p3.phi(), imass),
          p4Cartesian_(p3.x(), p3.y(), p3.z(), iEnergy),
          qx3_(integerCharge ? q * 3 : q),
          pdgId_(pdgId),
          status_(status) {}

    /// set internal cache
    inline void setCartesian() { p4Cartesian_ = p4Polar_; }

    /// electric charge
    int charge() const { return qx3_ / 3; }
    /// set electric charge
    void setCharge(Charge q) { qx3_ = q * 3; }
    /// electric charge
    int threeCharge() const { return qx3_; }
    /// set electric charge
    void setThreeCharge(Charge qx3) { qx3_ = qx3; }
    /// four-momentum Lorentz vector
    const LorentzVector& p4() const { return p4Cartesian_; }
    /// four-momentum Lorentz vector
    const PolarLorentzVector& polarP4() const { return p4Polar_; }
    /// spatial momentum vector
    Vector momentum() const { return p4Cartesian_.Vect(); }
    /// boost vector to boost a Lorentz vector
    /// to the particle center of mass system
    Vector boostToCM() const { return p4Cartesian_.BoostToCM(); }
    /// magnitude of momentum vector
    double p() const { return p4Cartesian_.P(); }
    /// energy
    double energy() const { return p4Cartesian_.E(); }
    /// transverse energy
    double et() const { return (pt() <= 0) ? 0 : p4Cartesian_.Et(); }
    /// transverse energy squared (use this for cuts)!
    double et2() const { return (pt() <= 0) ? 0 : p4Cartesian_.Et2(); }
    /// mass
    double mass() const { return p4Polar_.mass(); }
    /// mass squared
    double massSqr() const { return mass() * mass(); }
    /// transverse mass
    double mt() const { return p4Polar_.Mt(); }
    /// transverse mass squared
    double mtSqr() const { return p4Polar_.Mt2(); }
    /// x coordinate of momentum vector
    double px() const { return p4Cartesian_.Px(); }
    /// y coordinate of momentum vector
    double py() const { return p4Cartesian_.Py(); }
    /// z coordinate of momentum vector
    double pz() const { return p4Cartesian_.Pz(); }
    /// transverse momentum
    double pt() const { return p4Polar_.pt(); }
    /// momentum azimuthal angle
    double phi() const { return p4Polar_.phi(); }
    /// momentum polar angle
    double theta() const { return p4Cartesian_.Theta(); }
    /// momentum pseudorapidity
    double eta() const { return p4Polar_.eta(); }
    /// repidity
    double rapidity() const { return p4Polar_.Rapidity(); }
    /// repidity
    double y() const { return rapidity(); }
    /// set 4-momentum
    void setP4(const LorentzVector& p4) {
      p4Cartesian_ = p4;
      p4Polar_ = p4;
    }
    /// set 4-momentum
    void setP4(const PolarLorentzVector& p4) {
      p4Polar_ = p4;
      p4Cartesian_ = p4;
    }
    /// set particle mass
    void setMass(double m) {
      p4Polar_.SetM(m);
      setCartesian();
    }
    void setPz(double pz) {
      p4Cartesian_.SetPz(pz);
      p4Polar_ = p4Cartesian_;
    }
    /// vertex position
    const Point& vertex() const { return vertex_; }
    /// x coordinate of vertex position
    double vx() const { return vertex_.X(); }
    /// y coordinate of vertex position
    double vy() const { return vertex_.Y(); }
    /// z coordinate of vertex position
    double vz() const { return vertex_.Z(); }
    /// set vertex
    void setVertex(const Point& vertex) { vertex_ = vertex; }
    /// PDG identifier
    int pdgId() const { return pdgId_; }
    // set PDG identifier
    void setPdgId(int pdgId) { pdgId_ = pdgId; }
    /// status word
    int status() const { return status_; }
    /// set status word
    void setStatus(int status) { status_ = status; }
    /// set long lived flag
    void setLongLived() { status_ |= longLivedTag; }
    /// is long lived?
    bool longLived() const { return status_ & longLivedTag; }
    /// set mass constraint flag
    void setMassConstraint() { status_ |= massConstraintTag; }
    /// do mass constraint?
    bool massConstraint() const { return status_ & massConstraintTag; }

  private:
    static const unsigned int longLivedTag = 65536;
    static const unsigned int massConstraintTag = 131072;

  private:
    /// vertex position
    Point vertex_;

    /// four-momentum Lorentz vector
    PolarLorentzVector p4Polar_;
    /// internal cache for p4
    LorentzVector p4Cartesian_;

    /// electric charge
    Charge qx3_;

    /// PDG identifier
    int pdgId_;
    /// status word
    int status_;
  };

}  // namespace reco

#endif