File indexing completed on 2023-10-25 10:02:09
0001 #include "RecoTauTag/RecoTau/interface/PositionAtECalEntranceComputer.h"
0002
0003 #include "MagneticField/Engine/interface/MagneticField.h"
0004 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0005 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0006 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0007 #include "CommonTools/BaseParticlePropagator/interface/BaseParticlePropagator.h"
0008
0009 #include <cassert>
0010
0011 PositionAtECalEntranceComputer::PositionAtECalEntranceComputer(edm::ConsumesCollector&& cc, bool isPhase2)
0012 : bField_esToken_(cc.esConsumes<MagneticField, IdealMagneticFieldRecord>()),
0013 caloGeo_esToken_(cc.esConsumes<CaloGeometry, CaloGeometryRecord>()),
0014 bField_z_(-1.),
0015 isPhase2_(isPhase2) {}
0016
0017 PositionAtECalEntranceComputer::PositionAtECalEntranceComputer(edm::ConsumesCollector& cc, bool isPhase2)
0018 : bField_esToken_(cc.esConsumes<MagneticField, IdealMagneticFieldRecord>()),
0019 caloGeo_esToken_(cc.esConsumes<CaloGeometry, CaloGeometryRecord>()),
0020 bField_z_(-1.),
0021 isPhase2_(isPhase2) {}
0022
0023 PositionAtECalEntranceComputer::~PositionAtECalEntranceComputer() {}
0024
0025 void PositionAtECalEntranceComputer::beginEvent(const edm::EventSetup& es) {
0026 bField_z_ = es.getData(bField_esToken_).inTesla(GlobalPoint(0., 0., 0.)).z();
0027 if (isPhase2_) {
0028 recHitTools_.setGeometry(es.getData(caloGeo_esToken_));
0029 hgcalFace_z_ = recHitTools_.getPositionLayer(1).z();
0030 }
0031 }
0032
0033 reco::Candidate::Point PositionAtECalEntranceComputer::operator()(const reco::Candidate* particle,
0034 bool& success) const {
0035 assert(bField_z_ != -1.);
0036 reco::Candidate::Point position;
0037 BaseParticlePropagator propagator = BaseParticlePropagator(
0038 RawParticle(particle->p4(),
0039 math::XYZTLorentzVector(particle->vertex().x(), particle->vertex().y(), particle->vertex().z(), 0.),
0040 particle->charge()),
0041 0.,
0042 0.,
0043 bField_z_);
0044 if (!isPhase2_ || std::abs(particle->eta()) < ecalBarrelEndcapEtaBorder_) {
0045 propagator.propagateToEcalEntrance(false);
0046 } else {
0047 if (std::abs(particle->vertex().z()) >= hgcalFace_z_) {
0048 success = false;
0049 return position;
0050 }
0051 propagator.setPropagationConditions(152.6, hgcalFace_z_, false);
0052 propagator.propagate();
0053 }
0054 if (propagator.getSuccess() != 0) {
0055 position = propagator.particle().vertex().Vect();
0056 success = (std::abs(position.eta()) <= hgcalHfEtaBorder_);
0057 } else {
0058 success = false;
0059 }
0060 return position;
0061 }