CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATPackedCandidateProducer.cc
Go to the documentation of this file.
1 #include <string>
2 
3 
22 
23 /*#include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
24 #include "TrackingTools/GeomPropagators/interface/AnalyticalImpactPointExtrapolator.h"
25 #include "MagneticField/Engine/interface/MagneticField.h"
26 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
27 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
28 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
29 */
30 //#define CRAZYSORT
31 
32 namespace pat {
34  const static int qualityMap[8] = {1,0,1,1,4,4,5,6};
35 
37  public:
40 
41  virtual void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
42 
43  //sorting of cands to maximize the zlib compression
45  if (std::abs(i.charge()) == std::abs(j.charge())) {
46  if(i.charge()!=0){
47  if(i.pt() > minPtForTrackProperties_ and j.pt() <= minPtForTrackProperties_ ) return true;
48  if(i.pt() <= minPtForTrackProperties_ and j.pt() > minPtForTrackProperties_ ) return false;
49  }
50  if(i.vertexRef() == j.vertexRef())
51  return i.eta() > j.eta();
52  else
53  return i.vertexRef().key() < j.vertexRef().key();
54  }
55  return std::abs(i.charge()) > std::abs(j.charge());
56  }
57  template <typename T>
58  std::vector<size_t> sort_indexes(const std::vector<T> &v ) const {
59  std::vector<size_t> idx(v.size());
60  for (size_t i = 0; i != idx.size(); ++i) idx[i] = i;
61  std::sort(idx.begin(), idx.end(),[&v,this](size_t i1, size_t i2) { return candsOrdering(v[i1],v[i2]);});
62  return idx;
63  }
64 
65  private:
66  //if PuppiSrc && PuppiNoLepSrc are empty, usePuppi is false
67  //otherwise assumes that if they are set, you wanted to use puppi and will throw an exception
68  //if the puppis are not found
69  const bool usePuppi_;
70 
82  std::vector< edm::EDGetTokenT<edm::View<reco::CompositePtrCandidate> > > SVWhiteLists_;
83 
85 
86  // for debugging
87  float calcDxy(float dx, float dy, float phi) const {
88  return - dx * std::sin(phi) + dy * std::cos(phi);
89  }
91  return p.Z()-v.Z() - ((p.X()-v.X()) * c.px() + (p.Y()-v.Y())*c.py()) * c.pz()/(c.pt()*c.pt());
92  }
93  };
94 }
95 
97  usePuppi_(!iConfig.getParameter<edm::InputTag>("PuppiSrc").encode().empty() ||
98  !iConfig.getParameter<edm::InputTag>("PuppiNoLepSrc").encode().empty()),
99  Cands_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("inputCollection"))),
100  PVs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("inputVertices"))),
101  PVAsso_(consumes<edm::Association<reco::VertexCollection> >(iConfig.getParameter<edm::InputTag>("vertexAssociator"))),
102  PVAssoQuality_(consumes<edm::ValueMap<int> >(iConfig.getParameter<edm::InputTag>("vertexAssociator"))),
103  PVOrigs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("originalVertices"))),
104  TKOrigs_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("originalTracks"))),
105  PuppiWeight_(usePuppi_ ? consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("PuppiSrc")) : edm::EDGetTokenT< edm::ValueMap<float> >()),
106  PuppiWeightNoLep_(usePuppi_ ? consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("PuppiNoLepSrc")) : edm::EDGetTokenT< edm::ValueMap<float> >()),
107  PuppiCandsMap_(usePuppi_ ? consumes<edm::ValueMap<reco::CandidatePtr> >(iConfig.getParameter<edm::InputTag>("PuppiSrc")) : edm::EDGetTokenT<edm::ValueMap<reco::CandidatePtr> >() ),
108  PuppiCands_(usePuppi_ ? consumes<std::vector< reco::PFCandidate > >(iConfig.getParameter<edm::InputTag>("PuppiSrc")) : edm::EDGetTokenT<std::vector< reco::PFCandidate > >() ),
109  PuppiCandsNoLep_(usePuppi_ ? consumes<std::vector< reco::PFCandidate > >(iConfig.getParameter<edm::InputTag>("PuppiNoLepSrc")) : edm::EDGetTokenT<std::vector< reco::PFCandidate > >()),
110  minPtForTrackProperties_(iConfig.getParameter<double>("minPtForTrackProperties"))
111 {
112  std::vector<edm::InputTag> sv_tags = iConfig.getParameter<std::vector<edm::InputTag> >("secondaryVerticesForWhiteList");
113  for(auto itag : sv_tags){
114  SVWhiteLists_.push_back(
116  );
117  }
118 
119  produces< std::vector<pat::PackedCandidate> > ();
120  produces< edm::Association<pat::PackedCandidateCollection> > ();
121  produces< edm::Association<reco::PFCandidateCollection> > ();
122 }
123 
125 
126 
127 
129 
131  iEvent.getByToken( Cands_, cands );
132  std::vector<reco::Candidate>::const_iterator cand;
133 
134  edm::Handle<edm::ValueMap<float> > puppiWeight;
137  edm::Handle<edm::ValueMap<float> > puppiWeightNoLep;
139  std::vector<reco::CandidatePtr> puppiCandsNoLepPtrs;
140  if(usePuppi_){
141  iEvent.getByToken( PuppiWeight_, puppiWeight );
142  iEvent.getByToken( PuppiCandsMap_, puppiCandsMap );
143  iEvent.getByToken( PuppiCands_, puppiCands );
144  iEvent.getByToken( PuppiWeightNoLep_, puppiWeightNoLep );
145  iEvent.getByToken( PuppiCandsNoLep_, puppiCandsNoLep );
146  for (auto pup : *puppiCandsNoLep){
147  puppiCandsNoLepPtrs.push_back(pup.sourceCandidatePtr(0));
148  }
149  }
150  std::vector<int> mappingPuppi(usePuppi_ ? puppiCands->size() : 0);
151 
153  iEvent.getByToken( PVOrigs_, PVOrigs );
154 
156  iEvent.getByToken(PVAsso_,assoHandle);
157  edm::Handle<edm::ValueMap<int> > assoQualityHandle;
158  iEvent.getByToken(PVAssoQuality_,assoQualityHandle);
159  const edm::Association<reco::VertexCollection> & associatedPV=*(assoHandle.product());
160  const edm::ValueMap<int> & associationQuality=*(assoQualityHandle.product());
161 
162 
163  std::set<unsigned int> whiteList;
164  for(auto itoken : SVWhiteLists_) {
166  iEvent.getByToken(itoken, svWhiteListHandle);
167  const edm::View<reco::CompositePtrCandidate > & svWhiteList=*(svWhiteListHandle.product());
168  for(unsigned int i=0; i<svWhiteList.size();i++) {
169  for(unsigned int j=0; j< svWhiteList[i].numberOfSourceCandidatePtrs(); j++) {
170  const edm::Ptr<reco::Candidate> & c = svWhiteList[i].sourceCandidatePtr(j);
171  if(c.id() == cands.id()) whiteList.insert(c.key());
172  }
173  }
174  }
175 
176 
178  iEvent.getByToken( PVs_, PVs );
179  reco::VertexRef PV(PVs.id());
180  reco::VertexRefProd PVRefProd(PVs);
181  math::XYZPoint PVpos;
182 
183 
185  iEvent.getByToken( TKOrigs_, TKOrigs );
186  std::auto_ptr< std::vector<pat::PackedCandidate> > outPtrP( new std::vector<pat::PackedCandidate> );
187  std::vector<int> mapping(cands->size());
188  std::vector<int> mappingReverse(cands->size());
189  std::vector<int> mappingTk(TKOrigs->size(), -1);
190 
191  for(unsigned int ic=0, nc = cands->size(); ic < nc; ++ic) {
192  const reco::PFCandidate &cand=(*cands)[ic];
193  float phiAtVtx = cand.phi();
194  const reco::Track *ctrack = 0;
195  if ((abs(cand.pdgId()) == 11 || cand.pdgId() == 22) && cand.gsfTrackRef().isNonnull()) {
196  ctrack = &*cand.gsfTrackRef();
197  } else if (cand.trackRef().isNonnull()) {
198  ctrack = &*cand.trackRef();
199  }
200  if (ctrack) {
201  float dist=1e99;
202  int pvi=-1;
203  for(size_t ii=0;ii<PVs->size();ii++){
204  float dz=std::abs(ctrack->dz( ((*PVs)[ii]).position()));
205  if(dz<dist) {pvi=ii;dist=dz; }
206  }
207  PV = reco::VertexRef(PVs, pvi);
208  math::XYZPoint vtx = cand.vertex();
210  const reco::VertexRef & PVOrig = associatedPV[reco::CandidatePtr(cands,ic)];
211  if(PVOrig.isNonnull()) PV = reco::VertexRef(PVs, PVOrig.key()); // WARNING: assume the PV slimmer is keeping same order
212  int quality=associationQuality[reco::CandidatePtr(cands,ic)];
213 // if ((size_t)pvi!=PVOrig.key()) std::cout << "not closest in Z" << pvi << " " << PVOrig.key() << " " << cand.pt() << " " << quality << std::endl;
214  // TrajectoryStateOnSurface tsos = extrapolator.extrapolate(trajectoryStateTransform::initialFreeState(*ctrack,&*magneticField), RecoVertex::convertPos(PV->position()));
215  // vtx = tsos.globalPosition();
216  // phiAtVtx = tsos.globalDirection().phi();
217  vtx = ctrack->referencePoint();
218  phiAtVtx = ctrack->phi();
220  if (nlost == 0) {
221  if (ctrack->hitPattern().hasValidHitInFirstPixelBarrel()) {
223  }
224  } else {
226  }
227 
228 
229  outPtrP->push_back( pat::PackedCandidate(cand.polarP4(), vtx, phiAtVtx, cand.pdgId(), PVRefProd, PV.key()));
230  outPtrP->back().setAssociationQuality(pat::PackedCandidate::PVAssociationQuality(qualityMap[quality]));
231  if(cand.trackRef().isNonnull() && PVOrig->trackWeight(cand.trackRef()) > 0.5 && quality == 7) {
232  outPtrP->back().setAssociationQuality(pat::PackedCandidate::UsedInFitTight);
233  }
234  // properties of the best track
235  outPtrP->back().setLostInnerHits( lostHits );
236  if(outPtrP->back().pt() > minPtForTrackProperties_ || whiteList.find(ic)!=whiteList.end()) {
237  outPtrP->back().setTrackProperties(*ctrack);
238  //outPtrP->back().setTrackProperties(*ctrack,tsos.curvilinearError());
239  }
240 
241  // these things are always for the CKF track
242  outPtrP->back().setTrackHighPurity( cand.trackRef().isNonnull() && cand.trackRef()->quality(reco::Track::highPurity) );
243  if (cand.muonRef().isNonnull()) {
244  outPtrP->back().setMuonID(cand.muonRef()->isStandAloneMuon(), cand.muonRef()->isGlobalMuon());
245  }
246  } else {
247 
248  if (!PVs->empty()) {
249  PV = reco::VertexRef(PVs, 0);
250  PVpos = PV->position();
251  }
252 
253  outPtrP->push_back( pat::PackedCandidate(cand.polarP4(), PVpos, cand.phi(), cand.pdgId(), PVRefProd, PV.key()));
255  }
256 
257  // neutrals
258 
259  if(abs(cand.pdgId()) == 1 || abs(cand.pdgId()) == 130) {
260  outPtrP->back().setHcalFraction(cand.hcalEnergy()/(cand.ecalEnergy()+cand.hcalEnergy()));
261  } else {
262  outPtrP->back().setHcalFraction(0);
263  }
264 
265 
266  if (usePuppi_){
267  reco::PFCandidateRef pkref( cands, ic );
268  // outPtrP->back().setPuppiWeight( (*puppiWeight)[pkref]);
269 
270  float puppiWeightVal = (*puppiWeight)[pkref];
271  float puppiWeightNoLepVal = 0.0;
272  // Check the "no lepton" puppi weights.
273  // If present, then it is not a lepton, use stored weight
274  // If absent, it is a lepton, so set the weight to 1.0
275  if ( puppiWeightNoLep.isValid() ) {
276  // Look for the pointer inside the "no lepton" candidate collection.
277  auto pkrefPtr = pkref->sourceCandidatePtr(0);
278 
279  bool foundNoLep = false;
280  for ( size_t ipcnl = 0; ipcnl < puppiCandsNoLepPtrs.size(); ipcnl++){
281  if (puppiCandsNoLepPtrs[ipcnl] == pkrefPtr){
282  foundNoLep = true;
283  puppiWeightNoLepVal = puppiCandsNoLep->at(ipcnl).pt()/cand.pt(); // a hack for now, should use the value map
284  break;
285  }
286  }
287  if ( !foundNoLep || puppiWeightNoLepVal > 1 ) {
288  puppiWeightNoLepVal = 1.0;
289  }
290  }
291  outPtrP->back().setPuppiWeight( puppiWeightVal, puppiWeightNoLepVal );
292 
293  mappingPuppi[((*puppiCandsMap)[pkref]).key()]=ic;
294  }
295 
296  mapping[ic] = ic; // trivial at the moment!
297  if (cand.trackRef().isNonnull() && cand.trackRef().id() == TKOrigs.id()) {
298  mappingTk[cand.trackRef().key()] = ic;
299  }
300 
301  }
302 
303  std::auto_ptr< std::vector<pat::PackedCandidate> > outPtrPSorted( new std::vector<pat::PackedCandidate> );
304  std::vector<size_t> order=sort_indexes(*outPtrP);
305  std::vector<size_t> reverseOrder(order.size());
306  for(size_t i=0,nc=cands->size();i<nc;i++) {
307  outPtrPSorted->push_back((*outPtrP)[order[i]]);
308  reverseOrder[order[i]] = i;
309  mappingReverse[order[i]]=i;
310  }
311 
312  // Fix track association for sorted candidates
313  for(size_t i=0,ntk=mappingTk.size();i<ntk;i++){
314  if(mappingTk[i] >= 0)
315  mappingTk[i]=reverseOrder[mappingTk[i]];
316  }
317 
318  for(size_t i=0,ntk=mappingPuppi.size();i<ntk;i++){
319  mappingPuppi[i]=reverseOrder[mappingPuppi[i]];
320  }
321 
322  edm::OrphanHandle<pat::PackedCandidateCollection> oh = iEvent.put( outPtrPSorted );
323 
324  // now build the two maps
325  std::auto_ptr<edm::Association<pat::PackedCandidateCollection> > pf2pc(new edm::Association<pat::PackedCandidateCollection>(oh ));
326  std::auto_ptr<edm::Association<reco::PFCandidateCollection > > pc2pf(new edm::Association<reco::PFCandidateCollection >(cands));
329  pf2pcFiller.insert(cands, mappingReverse.begin(), mappingReverse.end());
330  pc2pfFiller.insert(oh , order.begin(), order.end());
331  // include also the mapping track -> packed PFCand
332  pf2pcFiller.insert(TKOrigs, mappingTk.begin(), mappingTk.end());
333  if(usePuppi_) pf2pcFiller.insert(puppiCands, mappingPuppi.begin(), mappingPuppi.end());
334 
335  pf2pcFiller.fill();
336  pc2pfFiller.fill();
337  iEvent.put(pf2pc);
338  iEvent.put(pc2pf);
339 
340 }
341 
342 
T getParameter(std::string const &) const
double ecalEnergy() const
return corrected Ecal energy
Definition: PFCandidate.h:221
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:676
int i
Definition: DBlmapReader.cc:9
const edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeight_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
key_type key() const
Definition: Ptr.h:186
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual double pz() const =0
z coordinate of momentum vector
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
virtual double phi() const final
momentum azimuthal angle
size_type size() const
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:640
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static const int qualityMap[8]
conversion map from quality flags used in PV association and miniAOD one
const edm::EDGetTokenT< reco::PFCandidateCollection > Cands_
virtual double pt() const
transverse momentum
key_type key() const
Accessor for product key.
Definition: Ref.h:264
int ii
Definition: cuy.py:588
virtual double eta() const
momentum pseudorapidity
const reco::VertexRef vertexRef() const
ProductID id() const
Accessor for product ID.
Definition: Ref.h:258
const edm::EDGetTokenT< edm::ValueMap< int > > PVAssoQuality_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
int iEvent
Definition: GenABIO.cc:230
PATPackedCandidateProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float calcDxy(float dx, float dy, float phi) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual int charge() const
electric charge
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:647
const edm::EDGetTokenT< reco::TrackCollection > TKOrigs_
const edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCands_
std::vector< size_t > sort_indexes(const std::vector< T > &v) const
virtual double py() const =0
y coordinate of momentum vector
const edm::EDGetTokenT< edm::ValueMap< reco::CandidatePtr > > PuppiCandsMap_
bool isValid() const
Definition: HandleBase.h:75
const edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeightNoLep_
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
reco::MuonRef muonRef() const
Definition: PFCandidate.cc:450
virtual double px() const =0
x coordinate of momentum vector
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:604
LostInnerHits
Enumerator specifying the.
std::vector< edm::EDGetTokenT< edm::View< reco::CompositePtrCandidate > > > SVWhiteLists_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
float calcDz(reco::Candidate::Point p, reco::Candidate::Point v, const reco::Candidate &c) const
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
T const * product() const
Definition: Handle.h:81
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
const edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCandsNoLep_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:181
const edm::EDGetTokenT< reco::VertexCollection > PVOrigs_
virtual const PolarLorentzVector & polarP4() const final
four-momentum Lorentz vector
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:902
virtual int pdgId() const final
PDG identifier.
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:471
double hcalEnergy() const
return corrected Hcal energy
Definition: PFCandidate.h:231
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
const edm::EDGetTokenT< edm::Association< reco::VertexCollection > > PVAsso_
bool hasValidHitInFirstPixelBarrel() const
Definition: HitPattern.cc:326
const edm::EDGetTokenT< reco::VertexCollection > PVs_
bool candsOrdering(pat::PackedCandidate i, pat::PackedCandidate j) const
virtual double pt() const final
transverse momentum