File indexing completed on 2024-04-06 12:04:51
0001 #ifndef DataFormats_ParticleFlowReco_PFCluster_h
0002 #define DataFormats_ParticleFlowReco_PFCluster_h
0003
0004 #include <iostream>
0005 #include <vector>
0006 #include <algorithm>
0007
0008 #include <Rtypes.h>
0009
0010 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0011 #include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
0012 #include "DataFormats/Math/interface/Point3D.h"
0013 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
0014 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0015 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
0016 #include "Math/GenVector/PositionVector3D.h"
0017
0018 class PFClusterAlgo;
0019
0020 namespace reco {
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042 class PFCluster : public CaloCluster {
0043 public:
0044 typedef std::vector<std::pair<CaloClusterPtr::key_type, edm::Ptr<PFCluster>>> EEtoPSAssociation;
0045
0046
0047
0048 typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double>> REPPoint;
0049
0050 PFCluster() : CaloCluster(CaloCluster::particleFlow), time_(-99.0), layer_(PFLayer::NONE) {}
0051
0052
0053 PFCluster(PFLayer::Layer layer, double energy, double x, double y, double z);
0054
0055
0056 void reset();
0057
0058
0059 void resetHitsAndFractions();
0060
0061
0062 void addRecHitFraction(const reco::PFRecHitFraction& frac);
0063
0064
0065 const std::vector<reco::PFRecHitFraction>& recHitFractions() const { return rechits_; }
0066
0067
0068 void setLayer(PFLayer::Layer layer);
0069
0070
0071 PFLayer::Layer layer() const;
0072
0073
0074 double energy() const { return energy_; }
0075
0076
0077 float time() const { return time_; }
0078
0079 float timeError() const { return timeError_; }
0080
0081
0082 double depth() const { return depth_; }
0083
0084 void setTime(float time, float timeError = 0) {
0085 time_ = time;
0086 timeError_ = timeError;
0087 }
0088 void setTimeError(float timeError) { timeError_ = timeError; }
0089 void setDepth(double depth) { depth_ = depth; }
0090
0091
0092 const REPPoint& positionREP() const { return posrep_; }
0093
0094
0095 void calculatePositionREP() { posrep_.SetCoordinates(position_.Rho(), position_.Eta(), position_.Phi()); }
0096
0097
0098 static double getDepthCorrection(double energy, bool isBelowPS = false, bool isHadron = false);
0099
0100
0101
0102
0103
0104
0105 double charge() const { return 0; }
0106
0107
0108 double pt() const { return (energy() * sin(position_.theta())); }
0109
0110
0111 double theta() const { return position_.theta(); }
0112
0113
0114 math::XYZPoint const& vertex() const { return dummyVtx_; }
0115 double vx() const { return vertex().x(); }
0116 double vy() const { return vertex().y(); }
0117 double vz() const { return vertex().z(); }
0118
0119 #if !defined(__CINT__) && !defined(__MAKECINT__) && !defined(__REFLEX__)
0120 template <typename pruner>
0121 void pruneUsing(pruner prune) {
0122
0123 auto iter = std::find_if_not(rechits_.begin(), rechits_.end(), prune);
0124 if (iter == rechits_.end())
0125 return;
0126 auto first = iter - rechits_.begin();
0127 for (auto i = first; ++i < int(rechits_.size());) {
0128 if (prune(rechits_[i])) {
0129 rechits_[first] = std::move(rechits_[i]);
0130 hitsAndFractions_[first] = std::move(hitsAndFractions_[i]);
0131 ++first;
0132 }
0133 }
0134 rechits_.erase(rechits_.begin() + first, rechits_.end());
0135 hitsAndFractions_.erase(hitsAndFractions_.begin() + first, hitsAndFractions_.end());
0136 }
0137 #endif
0138
0139 private:
0140
0141 std::vector<reco::PFRecHitFraction> rechits_;
0142
0143
0144 REPPoint posrep_;
0145
0146
0147 float time_;
0148 float timeError_{0};
0149 double depth_{0};
0150
0151
0152 PFLayer::Layer layer_;
0153
0154
0155 static const constexpr double depthCorA_ = 0.89;
0156 static const constexpr double depthCorB_ = 7.3;
0157 static const constexpr double depthCorAp_ = 0.89;
0158 static const constexpr double depthCorBp_ = 4.0;
0159
0160 static const math::XYZPoint dummyVtx_;
0161 };
0162
0163 std::ostream& operator<<(std::ostream& out, const PFCluster& cluster);
0164
0165 }
0166
0167 #endif