CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
pat::PATPackedCandidateProducer Class Reference
Inheritance diagram for pat::PATPackedCandidateProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

bool candsOrdering (pat::PackedCandidate i, pat::PackedCandidate j)
 
 PATPackedCandidateProducer (const edm::ParameterSet &)
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
template<typename T >
std::vector< size_t > sort_indexes (const std::vector< T > &v)
 
 ~PATPackedCandidateProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

float calcDxy (float dx, float dy, float phi)
 
float calcDz (reco::Candidate::Point p, reco::Candidate::Point v, const reco::Candidate &c)
 

Private Attributes

edm::EDGetTokenT
< reco::PFCandidateCollection
Cands_
 
double minPtForTrackProperties_
 
edm::EDGetTokenT< std::vector
< reco::PFCandidate > > 
PuppiCands_
 
edm::EDGetTokenT
< edm::ValueMap
< reco::CandidatePtr > > 
PuppiCandsMap_
 
edm::EDGetTokenT< std::vector
< reco::PFCandidate > > 
PuppiCandsNoLep_
 
edm::EDGetTokenT
< edm::ValueMap< float > > 
PuppiWeight_
 
edm::EDGetTokenT
< edm::ValueMap< float > > 
PuppiWeightNoLep_
 
edm::EDGetTokenT
< edm::Association
< reco::VertexCollection > > 
PVAsso_
 
edm::EDGetTokenT
< edm::ValueMap< int > > 
PVAssoQuality_
 
edm::EDGetTokenT
< reco::VertexCollection
PVOrigs_
 
edm::EDGetTokenT
< reco::VertexCollection
PVs_
 
edm::EDGetTokenT< edm::View
< reco::CompositePtrCandidate > > 
SVWhiteList_
 
edm::EDGetTokenT
< reco::TrackCollection
TKOrigs_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 36 of file PATPackedCandidateProducer.cc.

Constructor & Destructor Documentation

pat::PATPackedCandidateProducer::PATPackedCandidateProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 90 of file PATPackedCandidateProducer.cc.

90  :
91  Cands_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("inputCollection"))),
92  PVs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("inputVertices"))),
94  PVAssoQuality_(consumes<edm::ValueMap<int> >(iConfig.getParameter<edm::InputTag>("vertexAssociator"))),
95  PVOrigs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("originalVertices"))),
96  TKOrigs_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("originalTracks"))),
100  PuppiCands_(consumes<std::vector< reco::PFCandidate > >(iConfig.getParameter<edm::InputTag>("PuppiSrc"))),
101  PuppiCandsNoLep_(consumes<std::vector< reco::PFCandidate > >(iConfig.getParameter<edm::InputTag>("PuppiNoLepSrc"))),
102  SVWhiteList_(consumes<edm::View< reco::CompositePtrCandidate > >(iConfig.getParameter<edm::InputTag>("secondaryVerticesForWhiteList"))),
103  minPtForTrackProperties_(iConfig.getParameter<double>("minPtForTrackProperties"))
104 {
105  produces< std::vector<pat::PackedCandidate> > ();
106  produces< edm::Association<pat::PackedCandidateCollection> > ();
107  produces< edm::Association<reco::PFCandidateCollection> > ();
108 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeightNoLep_
edm::EDGetTokenT< edm::View< reco::CompositePtrCandidate > > SVWhiteList_
edm::EDGetTokenT< reco::TrackCollection > TKOrigs_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::ValueMap< reco::CandidatePtr > > PuppiCandsMap_
edm::EDGetTokenT< reco::PFCandidateCollection > Cands_
edm::EDGetTokenT< reco::VertexCollection > PVs_
edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCands_
edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCandsNoLep_
edm::EDGetTokenT< edm::ValueMap< int > > PVAssoQuality_
edm::EDGetTokenT< edm::Association< reco::VertexCollection > > PVAsso_
edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeight_
edm::EDGetTokenT< reco::VertexCollection > PVOrigs_
pat::PATPackedCandidateProducer::~PATPackedCandidateProducer ( )

Definition at line 110 of file PATPackedCandidateProducer.cc.

110 {}

Member Function Documentation

float pat::PATPackedCandidateProducer::calcDxy ( float  dx,
float  dy,
float  phi 
)
inlineprivate

Definition at line 81 of file PATPackedCandidateProducer.cc.

References funct::cos(), and funct::sin().

81  {
82  return - dx * std::sin(phi) + dy * std::cos(phi);
83  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Definition: DDAxes.h:10
float pat::PATPackedCandidateProducer::calcDz ( reco::Candidate::Point  p,
reco::Candidate::Point  v,
const reco::Candidate c 
)
inlineprivate

Definition at line 84 of file PATPackedCandidateProducer.cc.

References reco::Candidate::pt(), reco::Candidate::px(), reco::Candidate::py(), and reco::Candidate::pz().

84  {
85  return p.Z()-v.Z() - ((p.X()-v.X()) * c.px() + (p.Y()-v.Y())*c.py()) * c.pz()/(c.pt()*c.pt());
86  }
virtual double pt() const =0
transverse momentum
virtual double pz() const =0
z coordinate of momentum vector
virtual double py() const =0
y coordinate of momentum vector
virtual double px() const =0
x coordinate of momentum vector
bool pat::PATPackedCandidateProducer::candsOrdering ( pat::PackedCandidate  i,
pat::PackedCandidate  j 
)
inline

Definition at line 44 of file PATPackedCandidateProducer.cc.

References funct::abs(), pat::PackedCandidate::charge(), pat::PackedCandidate::eta(), edm::Ref< C, T, F >::key(), minPtForTrackProperties_, pat::PackedCandidate::pt(), and pat::PackedCandidate::vertexRef().

Referenced by sort_indexes().

44  {
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  }
virtual double pt() const
transverse momentum
key_type key() const
Accessor for product key.
Definition: Ref.h:266
virtual double eta() const
momentum pseudorapidity
const reco::VertexRef vertexRef() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int charge() const
electric charge
void pat::PATPackedCandidateProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDProducer.

Definition at line 114 of file PATPackedCandidateProducer.cc.

References funct::abs(), EnergyCorrector::c, reco::TrackBase::dz(), edm::helper::Filler< Map >::fill(), edm::Event::getByToken(), reco::PFCandidate::gsfTrackRef(), reco::HitPattern::hasValidHitInFirstPixelBarrel(), reco::TrackBase::highPurity, reco::TrackBase::hitPattern(), i, edm::HandleBase::id(), edm::Ptr< T >::id(), edm::Ref< C, T, F >::id(), cuy::ii, edm::helper::Filler< Map >::insert(), edm::Ref< C, T, F >::isNonnull(), edm::HandleBase::isValid(), j, edm::Ptr< T >::key(), edm::Ref< C, T, F >::key(), relval_steps::key, reco::HitPattern::MISSING_INNER_HITS, pat::PackedCandidate::moreLostInnerHits, reco::PFCandidate::muonRef(), pat::PackedCandidate::noLostInnerHits, reco::HitPattern::numberOfLostHits(), pat::PackedCandidate::oneLostInnerHit, reco::LeafCandidate::pdgId(), reco::LeafCandidate::phi(), reco::TrackBase::phi(), reco::LeafCandidate::polarP4(), edm::Handle< T >::product(), reco::LeafCandidate::pt(), edm::Event::put(), HLT_25ns14e33_v1_cff::quality, pat::qualityMap, reco::TrackBase::referencePoint(), edm::View< T >::size(), reco::PFCandidate::trackRef(), pat::PackedCandidate::UsedInFitTight, pat::PackedCandidate::validHitInFirstPixelBarrelLayer, and reco::PFCandidate::vertex().

114  {
115 
117  iEvent.getByToken( Cands_, cands );
118  std::vector<reco::Candidate>::const_iterator cand;
119 
120  edm::Handle< edm::ValueMap<float> > puppiWeight;
121  iEvent.getByToken( PuppiWeight_, puppiWeight );
123  iEvent.getByToken( PuppiCandsMap_, puppiCandsMap );
125  iEvent.getByToken( PuppiCands_, puppiCands );
126  std::vector<int> mappingPuppi(puppiCands->size());
127 
128  edm::Handle< edm::ValueMap<float> > puppiWeightNoLep;
129  iEvent.getByToken( PuppiWeightNoLep_, puppiWeightNoLep );
131  iEvent.getByToken( PuppiCandsNoLep_, puppiCandsNoLep );
132 
133  std::vector<reco::CandidatePtr> puppiCandsNoLepPtrs;
134  if (puppiCandsNoLep.isValid()){
135  for (auto pup : *puppiCandsNoLep){
136  puppiCandsNoLepPtrs.push_back(pup.sourceCandidatePtr(0));
137  }
138  }
139  auto const& puppiCandsNoLepV = puppiCandsNoLep.product();
140 
142  iEvent.getByToken( PVOrigs_, PVOrigs );
143 
145  iEvent.getByToken(PVAsso_,assoHandle);
146  edm::Handle<edm::ValueMap<int> > assoQualityHandle;
147  iEvent.getByToken(PVAssoQuality_,assoQualityHandle);
148  const edm::Association<reco::VertexCollection> & associatedPV=*(assoHandle.product());
149  const edm::ValueMap<int> & associationQuality=*(assoQualityHandle.product());
150 
152  iEvent.getByToken(SVWhiteList_,svWhiteListHandle);
153  const edm::View<reco::CompositePtrCandidate > & svWhiteList=*(svWhiteListHandle.product());
154  std::set<unsigned int> whiteList;
155  for(unsigned int i=0; i<svWhiteList.size();i++)
156  {
157  for(unsigned int j=0; j< svWhiteList[i].numberOfSourceCandidatePtrs(); j++) {
158  const edm::Ptr<reco::Candidate> & c = svWhiteList[i].sourceCandidatePtr(j);
159  if(c.id() == cands.id()) whiteList.insert(c.key());
160  }
161  }
162 
163 
165  iEvent.getByToken( PVs_, PVs );
166  reco::VertexRef PV(PVs.id());
167  math::XYZPoint PVpos;
168 
169 
171  iEvent.getByToken( TKOrigs_, TKOrigs );
172  std::auto_ptr< std::vector<pat::PackedCandidate> > outPtrP( new std::vector<pat::PackedCandidate> );
173  std::vector<int> mapping(cands->size());
174  std::vector<int> mappingReverse(cands->size());
175  std::vector<int> mappingTk(TKOrigs->size(), -1);
176 
177  for(unsigned int ic=0, nc = cands->size(); ic < nc; ++ic) {
178  const reco::PFCandidate &cand=(*cands)[ic];
179  float phiAtVtx = cand.phi();
180  const reco::Track *ctrack = 0;
181  if ((abs(cand.pdgId()) == 11 || cand.pdgId() == 22) && cand.gsfTrackRef().isNonnull()) {
182  ctrack = &*cand.gsfTrackRef();
183  } else if (cand.trackRef().isNonnull()) {
184  ctrack = &*cand.trackRef();
185  }
186  if (ctrack) {
187  float dist=1e99;
188  int pvi=-1;
189  for(size_t ii=0;ii<PVs->size();ii++){
190  float dz=std::abs(ctrack->dz( ((*PVs)[ii]).position()));
191  if(dz<dist) {pvi=ii;dist=dz; }
192  }
193  PV = reco::VertexRef(PVs, pvi);
194  math::XYZPoint vtx = cand.vertex();
196  const reco::VertexRef & PVOrig = associatedPV[reco::CandidatePtr(cands,ic)];
197  if(PVOrig.isNonnull()) PV = reco::VertexRef(PVs, PVOrig.key()); // WARNING: assume the PV slimmer is keeping same order
198  int quality=associationQuality[reco::CandidatePtr(cands,ic)];
199 // if ((size_t)pvi!=PVOrig.key()) std::cout << "not closest in Z" << pvi << " " << PVOrig.key() << " " << cand.pt() << " " << quality << std::endl;
200  // TrajectoryStateOnSurface tsos = extrapolator.extrapolate(trajectoryStateTransform::initialFreeState(*ctrack,&*magneticField), RecoVertex::convertPos(PV->position()));
201  // vtx = tsos.globalPosition();
202  // phiAtVtx = tsos.globalDirection().phi();
203  vtx = ctrack->referencePoint();
204  phiAtVtx = ctrack->phi();
206  if (nlost == 0) {
207  if (ctrack->hitPattern().hasValidHitInFirstPixelBarrel()) {
209  }
210  } else {
212  }
213 
214 
215  outPtrP->push_back( pat::PackedCandidate(cand.polarP4(), vtx, phiAtVtx, cand.pdgId(), PV));
216  outPtrP->back().setAssociationQuality(pat::PackedCandidate::PVAssociationQuality(qualityMap[quality]));
217  if(cand.trackRef().isNonnull() && PVOrig->trackWeight(cand.trackRef()) > 0.5 && quality == 7) {
218  outPtrP->back().setAssociationQuality(pat::PackedCandidate::UsedInFitTight);
219  }
220  // properties of the best track
221  outPtrP->back().setLostInnerHits( lostHits );
222  if(outPtrP->back().pt() > minPtForTrackProperties_ || whiteList.find(ic)!=whiteList.end()) {
223  outPtrP->back().setTrackProperties(*ctrack);
224  //outPtrP->back().setTrackProperties(*ctrack,tsos.curvilinearError());
225  }
226 
227  // these things are always for the CKF track
228  outPtrP->back().setTrackHighPurity( cand.trackRef().isNonnull() && cand.trackRef()->quality(reco::Track::highPurity) );
229  if (cand.muonRef().isNonnull()) {
230  outPtrP->back().setMuonID(cand.muonRef()->isStandAloneMuon(), cand.muonRef()->isGlobalMuon());
231  }
232  } else {
233 
234  if (!PVs->empty()) {
235  PV = reco::VertexRef(PVs, 0);
236  PVpos = PV->position();
237  }
238 
239  outPtrP->push_back( pat::PackedCandidate(cand.polarP4(), PVpos, cand.phi(), cand.pdgId(), PV));
241  }
242 
243  if (puppiWeight.isValid()){
244  reco::PFCandidateRef pkref( cands, ic );
245  // outPtrP->back().setPuppiWeight( (*puppiWeight)[pkref]);
246 
247  float puppiWeightVal = (*puppiWeight)[pkref];
248  float puppiWeightNoLepVal = 0.0;
249 
250  // Check the "no lepton" puppi weights.
251  // If present, then it is not a lepton, use stored weight
252  // If absent, it is a lepton, so set the weight to 1.0
253  if ( puppiWeightNoLep.isValid() ) {
254  // Look for the pointer inside the "no lepton" candidate collection.
255  auto pkrefPtr = pkref->sourceCandidatePtr(0);
256 
257  bool foundNoLep = false;
258  for ( size_t ipcnl = 0; ipcnl < puppiCandsNoLepPtrs.size(); ipcnl++){
259  if (puppiCandsNoLepPtrs[ipcnl] == pkrefPtr){
260  foundNoLep = true;
261  puppiWeightNoLepVal = puppiCandsNoLepV->at(ipcnl).pt()/cand.pt(); // a hack for now, should use the value map
262  break;
263  }
264  }
265  if ( !foundNoLep || puppiWeightNoLepVal > 1 ) {
266  puppiWeightNoLepVal = 1.0;
267  }
268  }
269  outPtrP->back().setPuppiWeight( puppiWeightVal, puppiWeightNoLepVal );
270 
271  mappingPuppi[((*puppiCandsMap)[pkref]).key()]=ic;
272  }
273 
274  mapping[ic] = ic; // trivial at the moment!
275  if (cand.trackRef().isNonnull() && cand.trackRef().id() == TKOrigs.id()) {
276  mappingTk[cand.trackRef().key()] = ic;
277  }
278 
279  }
280 
281  std::auto_ptr< std::vector<pat::PackedCandidate> > outPtrPSorted( new std::vector<pat::PackedCandidate> );
282  std::vector<size_t> order=sort_indexes(*outPtrP);
283  std::vector<size_t> reverseOrder(order.size());
284  for(size_t i=0,nc=cands->size();i<nc;i++) {
285  outPtrPSorted->push_back((*outPtrP)[order[i]]);
286  reverseOrder[order[i]] = i;
287  mappingReverse[order[i]]=i;
288  }
289 
290  // Fix track association for sorted candidates
291  for(size_t i=0,ntk=mappingTk.size();i<ntk;i++){
292  int origInd = mappingTk[i];
293  if (origInd>=0 ) mappingTk[i]=reverseOrder[origInd];
294  }
295 
296  for(size_t i=0,ntk=mappingPuppi.size();i<ntk;i++){
297  mappingPuppi[i]=reverseOrder[mappingPuppi[i]];
298  }
299 
300  edm::OrphanHandle<pat::PackedCandidateCollection> oh = iEvent.put( outPtrPSorted );
301 
302  // now build the two maps
303  std::auto_ptr<edm::Association<pat::PackedCandidateCollection> > pf2pc(new edm::Association<pat::PackedCandidateCollection>(oh ));
304  std::auto_ptr<edm::Association<reco::PFCandidateCollection > > pc2pf(new edm::Association<reco::PFCandidateCollection >(cands));
307  pf2pcFiller.insert(cands, mappingReverse.begin(), mappingReverse.end());
308  pc2pfFiller.insert(oh , order.begin(), order.end());
309  // include also the mapping track -> packed PFCand
310  pf2pcFiller.insert(TKOrigs, mappingTk.begin(), mappingTk.end());
311  pf2pcFiller.insert(puppiCands, mappingPuppi.begin(), mappingPuppi.end());
312 
313  pf2pcFiller.fill();
314  pc2pfFiller.fill();
315  iEvent.put(pf2pc);
316  iEvent.put(pc2pf);
317 
318 }
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:634
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeightNoLep_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
key_type key() const
Definition: Ptr.h:169
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
size_type size() const
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:598
key_type key() const
Accessor for product key.
Definition: Ref.h:266
int ii
Definition: cuy.py:588
edm::EDGetTokenT< edm::View< reco::CompositePtrCandidate > > SVWhiteList_
virtual double pt() const
transverse momentum
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
edm::EDGetTokenT< reco::TrackCollection > TKOrigs_
edm::EDGetTokenT< edm::ValueMap< reco::CandidatePtr > > PuppiCandsMap_
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
edm::EDGetTokenT< reco::PFCandidateCollection > Cands_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:647
edm::EDGetTokenT< reco::VertexCollection > PVs_
bool isValid() const
Definition: HandleBase.h:75
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCands_
reco::MuonRef muonRef() const
Definition: PFCandidate.cc:450
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:562
LostInnerHits
Enumerator specifying the.
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:411
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< size_t > sort_indexes(const std::vector< T > &v)
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:164
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:804
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
virtual const PolarLorentzVector & polarP4() const
four-momentum Lorentz vector
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:471
edm::EDGetTokenT< std::vector< reco::PFCandidate > > PuppiCandsNoLep_
static int qualityMap[8]
conversion map from quality flags used in PV association and miniAOD one
bool hasValidHitInFirstPixelBarrel() const
Definition: HitPattern.cc:276
edm::EDGetTokenT< edm::ValueMap< int > > PVAssoQuality_
edm::EDGetTokenT< edm::Association< reco::VertexCollection > > PVAsso_
edm::EDGetTokenT< edm::ValueMap< float > > PuppiWeight_
virtual double phi() const
momentum azimuthal angle
edm::EDGetTokenT< reco::VertexCollection > PVOrigs_
template<typename T >
std::vector<size_t> pat::PATPackedCandidateProducer::sort_indexes ( const std::vector< T > &  v)
inline

Definition at line 58 of file PATPackedCandidateProducer.cc.

References candsOrdering(), i, customizeTrackingMonitorSeedNumber::idx, python.multivaluedict::sort(), and findQualityFiles::v.

58  {
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  }
int i
Definition: DBlmapReader.cc:9
bool candsOrdering(pat::PackedCandidate i, pat::PackedCandidate j)
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...

Member Data Documentation

edm::EDGetTokenT<reco::PFCandidateCollection> pat::PATPackedCandidateProducer::Cands_
private

Definition at line 66 of file PATPackedCandidateProducer.cc.

double pat::PATPackedCandidateProducer::minPtForTrackProperties_
private

Definition at line 79 of file PATPackedCandidateProducer.cc.

Referenced by candsOrdering().

edm::EDGetTokenT<std::vector< reco::PFCandidate > > pat::PATPackedCandidateProducer::PuppiCands_
private

Definition at line 75 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<edm::ValueMap<reco::CandidatePtr> > pat::PATPackedCandidateProducer::PuppiCandsMap_
private

Definition at line 74 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<std::vector< reco::PFCandidate > > pat::PATPackedCandidateProducer::PuppiCandsNoLep_
private

Definition at line 76 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT< edm::ValueMap<float> > pat::PATPackedCandidateProducer::PuppiWeight_
private

Definition at line 72 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT< edm::ValueMap<float> > pat::PATPackedCandidateProducer::PuppiWeightNoLep_
private

Definition at line 73 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<edm::Association<reco::VertexCollection> > pat::PATPackedCandidateProducer::PVAsso_
private

Definition at line 68 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<edm::ValueMap<int> > pat::PATPackedCandidateProducer::PVAssoQuality_
private

Definition at line 69 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<reco::VertexCollection> pat::PATPackedCandidateProducer::PVOrigs_
private

Definition at line 70 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<reco::VertexCollection> pat::PATPackedCandidateProducer::PVs_
private

Definition at line 67 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<edm::View<reco::CompositePtrCandidate> > pat::PATPackedCandidateProducer::SVWhiteList_
private

Definition at line 77 of file PATPackedCandidateProducer.cc.

edm::EDGetTokenT<reco::TrackCollection> pat::PATPackedCandidateProducer::TKOrigs_
private

Definition at line 71 of file PATPackedCandidateProducer.cc.