CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelVertexCollectionTrimmer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoPixelVertexing/PixelVertexFinding
4 // Class: PixelVertexCollectionTrimmer
5 //
13 //
14 // Original Author: Riccardo Manzoni
15 // Created: Tue, 01 Apr 2014 10:11:16 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
36 
39 
41 public:
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46 
47 private:
48  virtual void produce(edm::Event&, const edm::EventSetup&) override;
49 
51  unsigned int maxVtx_ ;
52  double fractionSumPt2_ ;
53  double minSumPt2_ ;
54 
56 };
57 
59 {
60 
61  edm::InputTag vtxInputTag = iConfig.getParameter<edm::InputTag>("src" );
62  vtxToken_ = consumes<reco::VertexCollection>(vtxInputTag);
63  maxVtx_ = iConfig.getParameter<unsigned int> ("maxVtx" );
64  fractionSumPt2_ = iConfig.getParameter<double> ("fractionSumPt2" );
65  minSumPt2_ = iConfig.getParameter<double> ("minSumPt2" );
66 
67  edm::ParameterSet PVcomparerPSet = iConfig.getParameter<edm::ParameterSet>("PVcomparer");
68  double track_pt_min = PVcomparerPSet.getParameter<double>("track_pt_min");
69  double track_pt_max = PVcomparerPSet.getParameter<double>("track_pt_max");
70  double track_chi2_max = PVcomparerPSet.getParameter<double>("track_chi2_max");
71  double track_prob_min = PVcomparerPSet.getParameter<double>("track_prob_min");
72 
73  pvComparer_ = new PVClusterComparer(track_pt_min, track_pt_max, track_chi2_max, track_prob_min);
74 
75  produces<reco::VertexCollection>();
76 }
77 
79 {
80 }
81 
82 // ------------ method called to produce the data ------------
83 void
85 {
86  using namespace edm;
87 
88  std::auto_ptr<reco::VertexCollection> vtxs_trim(new reco::VertexCollection);
89 
91  iEvent.getByToken(vtxToken_,vtxs);
92 
93  double sumpt2 ;
94  //double sumpt2previous = -99. ;
95 
96  // this is not the logic we want, at least for now
97  // if requires the sumpt2 for vtx_n to be > threshold * sumpt2 vtx_n-1
98  // for (reco::VertexCollection::const_iterator vtx = vtxs->begin(); vtx != vtxs->end(); ++vtx, ++counter){
99  // if (counter > maxVtx_) break ;
100  // sumpt2 = PVCluster.pTSquaredSum(*vtx) ;
101  // if (sumpt2 > sumpt2previous*fractionSumPt2_ && sumpt2 > minSumPt2_ ) vtxs_trim->push_back(*vtx) ;
102  // else if (counter == 0 ) vtxs_trim->push_back(*vtx) ;
103  // sumpt2previous = sumpt2 ;
104  // }
105 
106  double sumpt2first = pvComparer_->pTSquaredSum(*(vtxs->begin())) ;
107 
108  for (reco::VertexCollection::const_iterator vtx = vtxs->begin(), evtx=vtxs->end();
109  vtx != evtx; ++vtx) {
110  if (vtxs_trim->size() >= maxVtx_) break ;
111  sumpt2 = pvComparer_->pTSquaredSum(*vtx) ;
112  // std::cout << "sumpt2: " << sumpt2 << "[" << sumpt2first << "]" << std::endl;
113  // if (sumpt2 >= sumpt2first*fractionSumPt2_ && sumpt2 > minSumPt2_ ) vtxs_trim->push_back(*vtx) ;
114  if (sumpt2 >= sumpt2first*fractionSumPt2_ && sumpt2 > minSumPt2_ ) vtxs_trim->push_back(*vtx) ;
115  }
116  // std::cout << " ==> # vertices: " << vtxs_trim->size() << std::endl;
117  iEvent.put(vtxs_trim);
118 
119 }
120 
121 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
122 void
124  //The following says we do not know what parameters are allowed so do no validation
125  // Please change this to state exactly what you do use, even if it is no parameters
127  desc.add<edm::InputTag> ("src", edm::InputTag(""))
128  ->setComment("input (pixel) vertex collection");
129  desc.add<unsigned int> ("maxVtx", 100)
130  ->setComment("max output collection size (number of accepted vertices)");
131  desc.add<double> ("fractionSumPt2", 0.3)
132  ->setComment("threshold on sumPt2 fraction of the leading vertex");
133  desc.add<double> ("minSumPt2", 0. )
134  ->setComment("min sumPt2");
135  edm::ParameterSetDescription PVcomparerPSet;
136  PVcomparerPSet.add<double>("track_pt_min",1.0)
137  ->setComment("min track p_T");
138  PVcomparerPSet.add<double>("track_pt_max",10.0)
139  ->setComment("max track p_T");
140  PVcomparerPSet.add<double>("track_chi2_max",99999.)
141  ->setComment("max track chi2");
142  PVcomparerPSet.add<double>("track_prob_min",-1.)
143  ->setComment("min track prob");
144  desc.add<edm::ParameterSetDescription>("PVcomparer", PVcomparerPSet)
145  ->setComment("from RecoPixelVertexing/PixelVertexFinding/python/PVClusterComparer_cfi.py");
146  descriptions.add("hltPixelVertexCollectionTrimmer",desc);
147 }
148 
149 //define this as a plug-in
T getParameter(std::string const &) const
void setComment(std::string const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
double pTSquaredSum(const PVCluster &v)
Calculate sum of square of the pT&#39;s of the tracks in the vertex.
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
virtual void produce(edm::Event &, const edm::EventSetup &) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
PixelVertexCollectionTrimmer(const edm::ParameterSet &)