CMS 3D CMS Logo

IsolatedPixelTrackCandidateProducer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <memory>
3 #include <algorithm>
4 #include <cmath>
5 
6 // Class header file
9 // Framework
17 //
18 
19 
20 // Math
21 #include "Math/GenVector/VectorUtil.h"
22 #include "Math/GenVector/PxPyPzE4D.h"
24 
28 
29 //magF
33 
34 //for ECAL geometry
40 
41 //#define DebugLog
42 
44  tok_hlt_( consumes<trigger::TriggerFilterObjectWithRefs>(config.getParameter<edm::InputTag>("L1GTSeedLabel")) ),
45  tok_l1_( consumes<l1extra::L1JetParticleCollection>(config.getParameter<edm::InputTag>("L1eTauJetsSource")) ),
46  tok_vert_( consumes<reco::VertexCollection>(config.getParameter<edm::InputTag>("VertexLabel")) ),
47  toks_pix_( edm::vector_transform(
48  config.getParameter<std::vector<edm::InputTag> >("PixelTracksSources"),
49  [this](edm::InputTag const & tag){return consumes<reco::TrackCollection>(tag);}) ),
50  bfield_( config.getParameter<std::string>("MagFieldRecordName") ),
51  prelimCone_( config.getParameter<double>("ExtrapolationConeSize") ),
52  pixelIsolationConeSizeAtEC_( config.getParameter<double>("PixelIsolationConeSizeAtEC") ),
53  vtxCutSeed_( config.getParameter<double>("MaxVtxDXYSeed") ),
54  vtxCutIsol_( config.getParameter<double>("MaxVtxDXYIsol") ),
55  tauAssocCone_( config.getParameter<double>("tauAssociationCone") ),
56  tauUnbiasCone_( config.getParameter<double>("tauUnbiasCone") ),
57  minPTrackValue_( config.getParameter<double>("minPTrack") ),
58  maxPForIsolationValue_( config.getParameter<double>("maxPTrackForIsolation") ),
59  ebEtaBoundary_( config.getParameter<double>("EBEtaBoundary") ),
60  rEB_( -1 ),
61  zEE_( -1 ),
62  bfVal_( 0 )
63 {
64  // Register the product
65  produces< reco::IsolatedPixelTrackCandidateCollection >();
66 }
67 
69 
72  std::vector<edm::InputTag> tracksrc = {edm::InputTag("hltPixelTracks")};
73  desc.add<edm::InputTag>("L1eTauJetsSource",edm::InputTag("hltCaloStage2Digis","Tau"));
74  desc.add<double>("tauAssociationCone", 0.0 );
75  desc.add<double>("tauUnbiasCone", 1.2 );
76  desc.add<std::vector<edm::InputTag> >("PixelTracksSources",tracksrc);
77  desc.add<double>("ExtrapolationConeSize", 1.0);
78  desc.add<double>("PixelIsolationConeSizeAtEC",40);
79  desc.add<edm::InputTag>("L1GTSeedLabel",edm::InputTag("hltL1sIsoTrack"));
80  desc.add<double>("MaxVtxDXYSeed", 101.0);
81  desc.add<double>("MaxVtxDXYIsol", 101.0);
82  desc.add<edm::InputTag>("VertexLabel",edm::InputTag("hltTrimmedPixelVertices"));
83  desc.add<std::string>("MagFieldRecordName","VolumeBasedMagneticField");
84  desc.add<double>("minPTrack", 5.0);
85  desc.add<double>("maxPTrackForIsolation", 3.0);
86  desc.add<double>("EBEtaBoundary", 1.479);
87  descriptions.add("isolPixelTrackProd",desc);
88 }
89 
91 
93  theEventSetup.get<CaloGeometryRecord>().get(pG);
94 
95  const double rad (dynamic_cast<const EcalBarrelGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalBarrel ))->avgRadiusXYFrontFaceCenter() ) ;
96  const double zz (dynamic_cast<const EcalEndcapGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalEndcap ))->avgAbsZFrontFaceCenter() ) ;
97 
98  rEB_ = rad;
99  zEE_ = zz;
100 
102  theEventSetup.get<IdealMagneticFieldRecord>().get(vbfField);
103  const VolumeBasedMagneticField* vbfCPtr = dynamic_cast<const VolumeBasedMagneticField*>(&(*vbfField));
104  GlobalVector BField = vbfCPtr->inTesla(GlobalPoint(0,0,0));
105  bfVal_=BField.mag();
106 }
107 
109 
110  auto trackCollection = std::make_unique<reco::IsolatedPixelTrackCandidateCollection>();
111 
112  //create vector of refs from input collections
113  std::vector<reco::TrackRef> pixelTrackRefs;
114 #ifdef DebugLog
115  edm::LogInfo("HcalIsoTrack") << "IsolatedPixelTrakCandidate: with" << toks_pix_.size() << " candidates to start with\n";
116 #endif
117  for (unsigned int iPix=0; iPix<toks_pix_.size(); iPix++) {
119  theEvent.getByToken(toks_pix_[iPix],iPixCol);
120  for (reco::TrackCollection::const_iterator pit=iPixCol->begin(); pit!=iPixCol->end(); pit++) {
121  pixelTrackRefs.push_back(reco::TrackRef(iPixCol,pit-iPixCol->begin()));
122  }
123  }
124 
126  theEvent.getByToken(tok_l1_,l1eTauJets);
127 
129  theEvent.getByToken(tok_vert_,pVert);
130 
131  double drMaxL1Track_ = tauAssocCone_;
132 
133  int ntr = 0;
134  std::vector<seedAtEC> VecSeedsatEC;
135  //loop to select isolated tracks
136  for (unsigned iS=0; iS<pixelTrackRefs.size(); iS++) {
137  bool vtxMatch = false;
138  //associate to vertex (in Z)
139  reco::VertexCollection::const_iterator vitSel;
140  double minDZ = 1000;
141  bool found(false);
142  for (reco::VertexCollection::const_iterator vit=pVert->begin(); vit!=pVert->end(); vit++) {
143  if (std::abs(pixelTrackRefs[iS]->dz(vit->position()))<minDZ) {
144  minDZ = std::abs(pixelTrackRefs[iS]->dz(vit->position()));
145  vitSel = vit;
146  found = true;
147  }
148  }
149  //cut on dYX:
150  if (found) {
151  if(std::abs(pixelTrackRefs[iS]->dxy(vitSel->position()))<vtxCutSeed_) vtxMatch=true;
152  } else {
153  vtxMatch=true;
154  }
155 
156  //check taujet matching
157  bool tmatch=false;
158  l1extra::L1JetParticleCollection::const_iterator selj;
159  for (l1extra::L1JetParticleCollection::const_iterator tj=l1eTauJets->begin(); tj!=l1eTauJets->end(); tj++) {
160  if (reco::deltaR(pixelTrackRefs[iS]->momentum().eta(), pixelTrackRefs[iS]->momentum().phi(), tj->momentum().eta(), tj->momentum().phi()) > drMaxL1Track_) continue;
161  selj = tj;
162  tmatch = true;
163  } //loop over L1 tau
164 
165  //propagate seed track to ECAL surface:
166  std::pair<double,double> seedCooAtEC;
167  // in case vertex is found:
168  if (found) seedCooAtEC=GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(), pixelTrackRefs[iS]->phi(), pixelTrackRefs[iS]->pt(), pixelTrackRefs[iS]->charge(), vitSel->z());
169  //in case vertex is not found:
170  else seedCooAtEC=GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(), pixelTrackRefs[iS]->phi(), pixelTrackRefs[iS]->pt(), pixelTrackRefs[iS]->charge(), 0);
171  seedAtEC seed(iS,(tmatch||vtxMatch),seedCooAtEC.first,seedCooAtEC.second);
172  VecSeedsatEC.push_back(seed);
173  }
174 #ifdef DebugLog
175  edm::LogInfo("HcalIsoTrack") << "IsolatedPixelTrakCandidate: " << VecSeedsatEC.size() << " seeds after propagation\n";
176 #endif
177 
178  for (unsigned int i=0; i<VecSeedsatEC.size(); i++) {
179  unsigned int iSeed = VecSeedsatEC[i].index;
180  if (!VecSeedsatEC[i].ok) continue;
181  if(pixelTrackRefs[iSeed]->p()<minPTrackValue_) continue;
182  l1extra::L1JetParticleCollection::const_iterator selj;
183  for (l1extra::L1JetParticleCollection::const_iterator tj=l1eTauJets->begin(); tj!=l1eTauJets->end(); tj++) {
184  if (reco::deltaR(pixelTrackRefs[iSeed]->momentum().eta(),pixelTrackRefs[iSeed]->momentum().phi(),tj->momentum().eta(),tj->momentum().phi()) > drMaxL1Track_) continue;
185  selj = tj;
186  } //loop over L1 tau
187  double maxP = 0;
188  double sumP = 0;
189  for(unsigned int j=0; j<VecSeedsatEC.size(); j++) {
190  if (i==j) continue;
191  unsigned int iSurr = VecSeedsatEC[j].index;
192  //define preliminary cone around seed track impact point from which tracks will be extrapolated:
193  if (reco::deltaR(pixelTrackRefs[iSeed]->eta(), pixelTrackRefs[iSeed]->phi(), pixelTrackRefs[iSurr]->eta(), pixelTrackRefs[iSurr]->phi())>prelimCone_) continue;
194  double minDZ2(1000);
195  bool found(false);
196  reco::VertexCollection::const_iterator vitSel2;
197  for (reco::VertexCollection::const_iterator vit=pVert->begin(); vit!=pVert->end(); vit++) {
198  if (std::abs(pixelTrackRefs[iSurr]->dz(vit->position()))<minDZ2) {
199  minDZ2 = std::abs(pixelTrackRefs[iSurr]->dz(vit->position()));
200  vitSel2 = vit;
201  found = true;
202  }
203  }
204  //cut ot dXY:
205  if (found&&std::abs(pixelTrackRefs[iSurr]->dxy(vitSel2->position()))>vtxCutIsol_) continue;
206  //calculate distance at ECAL surface and update isolation:
207  if (getDistInCM(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi, VecSeedsatEC[j].eta, VecSeedsatEC[j].phi)<pixelIsolationConeSizeAtEC_) {
208  sumP+=pixelTrackRefs[iSurr]->p();
209  if(pixelTrackRefs[iSurr]->p()>maxP) maxP=pixelTrackRefs[iSurr]->p();
210  }
211  }
212  if (maxP<maxPForIsolationValue_) {
213  reco::IsolatedPixelTrackCandidate newCandidate(pixelTrackRefs[iSeed], l1extra::L1JetParticleRef(l1eTauJets,selj-l1eTauJets->begin()), maxP, sumP);
214  newCandidate.setEtaPhiEcal(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi);
215  trackCollection->push_back(newCandidate);
216  ntr++;
217  }
218  }
219  // put the product in the event
220  theEvent.put(std::move(trackCollection));
221 #ifdef DebugLog
222  edm::LogInfo("HcalIsoTrack") << "IsolatedPixelTrackCandidate: Final # of candiates " << ntr << "\n";
223 #endif
224 }
225 
226 
227 double IsolatedPixelTrackCandidateProducer::getDistInCM(double eta1, double phi1, double eta2, double phi2) {
228  double Rec;
229  double theta1=2*atan(exp(-eta1));
230  double theta2=2*atan(exp(-eta2));
231  if (std::abs(eta1)<1.479) Rec=rEB_; //radius of ECAL barrel
232  else if (std::abs(eta1)>1.479&&std::abs(eta1)<7.0) Rec=tan(theta1)*zEE_; //distance from IP to ECAL endcap
233  else return 1000;
234 
235  //|vect| times tg of acos(scalar product)
236  double angle=acos((sin(theta1)*sin(theta2)*(sin(phi1)*sin(phi2)+cos(phi1)*cos(phi2))+cos(theta1)*cos(theta2)));
237  if (angle<M_PI_2) return std::abs((Rec/sin(theta1))*tan(angle));
238  else return 1000;
239 }
240 
241 
242 std::pair<double,double> IsolatedPixelTrackCandidateProducer::GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ) {
243 
244  double deltaPhi=0;
245  double etaEC = 100;
246  double phiEC = 100;
247 
248  double Rcurv = 9999999;
249  if (bfVal_!=0) Rcurv=pT*33.3*100/(bfVal_*10); //r(m)=pT(GeV)*33.3/B(kG)
250 
251  double ecDist = zEE_; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
252  double ecRad = rEB_; //radius of ECAL barrel (cm)
253  double theta = 2*atan(exp(-etaIP));
254  double zNew = 0;
255  if (theta>M_PI_2) theta=M_PI-theta;
256  if (std::abs(etaIP)<ebEtaBoundary_) {
257  if ((0.5*ecRad/Rcurv)>1) {
258  etaEC = 10000;
259  deltaPhi = 0;
260  } else {
261  deltaPhi =-charge*asin(0.5*ecRad/Rcurv);
262  double alpha1 = 2*asin(0.5*ecRad/Rcurv);
263  double z = ecRad/tan(theta);
264  if (etaIP>0) zNew = z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
265  else zNew =-z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
266  double zAbs=std::abs(zNew);
267  if (zAbs<ecDist) {
268  etaEC = -log(tan(0.5*atan(ecRad/zAbs)));
269  deltaPhi = -charge*asin(0.5*ecRad/Rcurv);
270  }
271  if (zAbs>ecDist) {
272  zAbs = (std::abs(etaIP)/etaIP)*ecDist;
273  double Zflight = std::abs(zAbs-vtxZ);
274  double alpha = (Zflight*ecRad)/(z*Rcurv);
275  double Rec = 2*Rcurv*sin(alpha/2);
276  deltaPhi =-charge*alpha/2;
277  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
278  }
279  }
280  } else {
281  zNew = (std::abs(etaIP)/etaIP)*ecDist;
282  double Zflight = std::abs(zNew-vtxZ);
283  double Rvirt = std::abs(Zflight*tan(theta));
284  double Rec = 2*Rcurv*sin(Rvirt/(2*Rcurv));
285  deltaPhi =-(charge)*(Rvirt/(2*Rcurv));
286  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
287  }
288 
289  if (zNew<0) etaEC=-etaEC;
290  phiEC = phiIP+deltaPhi;
291 
292  if (phiEC<-M_PI) phiEC += M_2_PI;
293  if (phiEC>M_PI) phiEC -= M_2_PI;
294 
295  std::pair<double,double> retVal(etaEC,phiEC);
296  return retVal;
297 }
298 
299 
std::pair< double, double > GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
float alpha
Definition: AMPTWrapper.h:95
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
void beginRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
GlobalVector inTesla(const GlobalPoint &g) const override
Field value ad specified global point, in Tesla.
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
IsolatedPixelTrackCandidateProducer(const edm::ParameterSet &ps)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
#define M_PI_2
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Definition: config.py:1
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T mag() const
Definition: PV3DBase.h:67
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setEtaPhiEcal(double eta, double phi)
eta, phi at ECAL surface
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< reco::VertexCollection > tok_vert_
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:63
double getDistInCM(double eta1, double phi1, double eta2, double phi2)
void produce(edm::Event &evt, const edm::EventSetup &es) override
def move(src, dest)
Definition: eostools.py:510
Definition: Run.h:44
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11