CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATVertexSlimmer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <memory>
3 
4 // user include files
9 
15 
16 namespace pat {
18  public:
19  explicit PATVertexSlimmer(const edm::ParameterSet&);
21 
22  virtual void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const;
23  private:
26  const bool rekeyScores_;
27  };
28 }
29 
31  src_(consumes<std::vector<reco::Vertex> >(iConfig.getParameter<edm::InputTag>("src"))),
32  score_(mayConsume<edm::ValueMap<float>>(iConfig.existsAs<edm::InputTag>("score")?iConfig.getParameter<edm::InputTag>("score"):edm::InputTag())),
33  rekeyScores_(iConfig.existsAs<edm::InputTag>("score"))
34 {
35  produces<std::vector<reco::Vertex> >();
36  if(rekeyScores_) produces<edm::ValueMap<float> >();
37 }
38 
40 
43  iEvent.getByToken(src_, vertices);
44  std::auto_ptr<std::vector<reco::Vertex> > outPtr(new std::vector<reco::Vertex>());
45 
46  outPtr->reserve(vertices->size());
47  for (unsigned int i = 0, n = vertices->size(); i < n; ++i) {
48  const reco::Vertex &v = (*vertices)[i];
49  auto co = v.covariance();
50  if(i>0) {
51  for(size_t j=0;j<3;j++){
52  for(size_t k=j;k<3;k++){
53  co(j,k) = MiniFloatConverter::reduceMantissaToNbits<10>( co(j,k) );
54  }
55  }
56  }
57  outPtr->push_back(reco::Vertex(v.position(), co, v.chi2(), v.ndof(), 0));
58  }
59 
60  auto oh = iEvent.put(outPtr);
61  if(rekeyScores_) {
63  iEvent.getByToken(score_, scores);
64  std::auto_ptr<edm::ValueMap<float> > vertexScoreOutput( new edm::ValueMap<float> );
65  edm::ValueMap<float>::const_iterator idIt=scores->begin();
66  for(;idIt!=scores->end();idIt++) {
67  if(idIt.id() == vertices.id()) break;
68  }
69  // std::find_if(scores->begin(), scores->end(), [vertices] (const edm::ValueMap<float>::const_iterator& s) { return s.id() == vertices.id(); } );
70  edm::ValueMap<float>::Filler vertexScoreFiller(*vertexScoreOutput);
71  vertexScoreFiller.insert(oh,idIt.begin(),idIt.end());
72  vertexScoreFiller.fill();
73  iEvent.put( vertexScoreOutput );
74  }
75 }
76 
int i
Definition: DBlmapReader.cc:9
PATVertexSlimmer(const edm::ParameterSet &)
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:123
const edm::EDGetTokenT< std::vector< reco::Vertex > > src_
const Point & position() const
position
Definition: Vertex.h:106
int iEvent
Definition: GenABIO.cc:230
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
double chi2() const
chi-squares
Definition: Vertex.h:95
int j
Definition: DBlmapReader.cc:9
double ndof() const
Definition: Vertex.h:102
const edm::EDGetTokenT< edm::ValueMap< float > > score_
container::const_iterator begin() const
Definition: ValueMap.h:176
ProductID id() const
Definition: ValueMap.h:175
container::const_iterator end() const
Definition: ValueMap.h:179