File indexing completed on 2024-04-06 12:27:27
0001 #include "RecoParticleFlow/PFProducer/interface/BlockElementImporterBase.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementBrem.h"
0005 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrack.h"
0006 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0007 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0008 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
0009 #include "RecoParticleFlow/PFProducer/interface/PFBlockElementSCEqual.h"
0010
0011 class GSFTrackImporter : public BlockElementImporterBase {
0012 public:
0013 GSFTrackImporter(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0014 : BlockElementImporterBase(conf, cc),
0015 _src(cc.consumes<reco::GsfPFRecTrackCollection>(conf.getParameter<edm::InputTag>("source"))),
0016 _isSecondary(conf.getParameter<bool>("gsfsAreSecondary")),
0017 _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {}
0018
0019 void importToBlock(const edm::Event&, ElementList&) const override;
0020
0021 private:
0022 edm::EDGetTokenT<reco::GsfPFRecTrackCollection> _src;
0023 const bool _isSecondary, _superClustersArePF;
0024 };
0025
0026 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, GSFTrackImporter, "GSFTrackImporter");
0027
0028 void GSFTrackImporter::importToBlock(const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
0029 auto gsftracks = e.getHandle(_src);
0030 elems.reserve(elems.size() + gsftracks->size());
0031
0032 auto SCs_end =
0033 std::partition(elems.begin(), elems.end(), [](const auto& a) { return a->type() == reco::PFBlockElement::SC; });
0034
0035 auto bgsf = gsftracks->cbegin();
0036 auto egsf = gsftracks->cend();
0037 for (auto gsftrack = bgsf; gsftrack != egsf; ++gsftrack) {
0038 reco::GsfPFRecTrackRef gsfref(gsftracks, std::distance(bgsf, gsftrack));
0039 const reco::GsfTrackRef& basegsfref = gsftrack->gsfTrackRef();
0040 const auto& gsfextraref = basegsfref->extra();
0041
0042 if (gsfextraref.isAvailable() && gsfextraref->seedRef().isAvailable()) {
0043 reco::ElectronSeedRef seedref = gsfextraref->seedRef().castTo<reco::ElectronSeedRef>();
0044 if (seedref.isAvailable() && seedref->isEcalDriven()) {
0045 reco::SuperClusterRef scref = seedref->caloCluster().castTo<reco::SuperClusterRef>();
0046 if (scref.isNonnull()) {
0047
0048 if (scref->seed()->seed().det() == DetId::HGCalEE || scref->seed()->seed().det() == DetId::HGCalHSi ||
0049 scref->seed()->seed().det() == DetId::HGCalHSc || scref->seed()->seed().det() == DetId::Forward)
0050 continue;
0051 PFBlockElementSCEqual myEqual(scref);
0052 auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
0053 if (sc_elem != SCs_end) {
0054 reco::PFBlockElementSuperCluster* scbe = static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
0055 scbe->setFromGsfElectron(true);
0056 } else {
0057 reco::PFBlockElementSuperCluster* scbe = new reco::PFBlockElementSuperCluster(scref);
0058 scbe->setFromGsfElectron(true);
0059 scbe->setFromPFSuperCluster(_superClustersArePF);
0060 SCs_end = elems.emplace(SCs_end, scbe);
0061 ++SCs_end;
0062 }
0063 }
0064 }
0065 }
0066
0067 size_t SCs_end_position = std::distance(elems.begin(), SCs_end);
0068
0069 const std::vector<reco::PFTrajectoryPoint>& PfGsfPoint = gsftrack->trajectoryPoints();
0070 unsigned int c_gsf = 0;
0071 bool PassTracker = false;
0072 unsigned int IndexPout = 0;
0073 for (const auto& pfGsfPoint : PfGsfPoint) {
0074 if (pfGsfPoint.isValid()) {
0075 int layGsfP = pfGsfPoint.layer();
0076 if (layGsfP == -1)
0077 PassTracker = true;
0078 else if (PassTracker && layGsfP > 0) {
0079 IndexPout = c_gsf - 1;
0080 break;
0081 }
0082 ++c_gsf;
0083 }
0084 }
0085 const math::XYZTLorentzVector& pin = PfGsfPoint[0].momentum();
0086 const math::XYZTLorentzVector& pout = PfGsfPoint[IndexPout].momentum();
0087 reco::PFBlockElementGsfTrack* temp = new reco::PFBlockElementGsfTrack(gsfref, pin, pout);
0088 if (_isSecondary) {
0089 temp->setTrackType(reco::PFBlockElement::T_FROM_GAMMACONV, true);
0090 }
0091 elems.emplace_back(temp);
0092
0093 for (const auto& brem : gsfref->PFRecBrem()) {
0094 const unsigned TrajP = brem.indTrajPoint();
0095 if (TrajP != 99) {
0096 const double DP = brem.DeltaP();
0097 const double sDP = brem.SigmaDeltaP();
0098 elems.emplace_back(new reco::PFBlockElementBrem(gsfref, DP, sDP, TrajP));
0099 }
0100 }
0101
0102 SCs_end = elems.begin() + SCs_end_position;
0103 }
0104 elems.shrink_to_fit();
0105 }