CMS 3D CMS Logo

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 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
30 
35 
38 
40 public:
43 
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46 private:
47  void produce(edm::Event&, const edm::EventSetup&) override;
48 
50  unsigned int maxVtx_;
52  double minSumPt2_;
53 
55 };
56 
59  vtxToken_ = consumes<reco::VertexCollection>(vtxInputTag);
60  maxVtx_ = iConfig.getParameter<unsigned int>("maxVtx");
61  fractionSumPt2_ = iConfig.getParameter<double>("fractionSumPt2");
62  minSumPt2_ = iConfig.getParameter<double>("minSumPt2");
63 
64  edm::ParameterSet PVcomparerPSet = iConfig.getParameter<edm::ParameterSet>("PVcomparer");
65  double track_pt_min = PVcomparerPSet.getParameter<double>("track_pt_min");
66  double track_pt_max = PVcomparerPSet.getParameter<double>("track_pt_max");
67  double track_chi2_max = PVcomparerPSet.getParameter<double>("track_chi2_max");
68  double track_prob_min = PVcomparerPSet.getParameter<double>("track_prob_min");
69 
71 
72  produces<reco::VertexCollection>();
73 }
74 
76 
77 // ------------ method called to produce the data ------------
79  using namespace edm;
80 
81  auto vtxs_trim = std::make_unique<reco::VertexCollection>();
82 
84  iEvent.getByToken(vtxToken_, vtxs);
85 
86  double sumpt2;
87  //double sumpt2previous = -99. ;
88 
89  // this is not the logic we want, at least for now
90  // if requires the sumpt2 for vtx_n to be > threshold * sumpt2 vtx_n-1
91  // for (reco::VertexCollection::const_iterator vtx = vtxs->begin(); vtx != vtxs->end(); ++vtx, ++counter){
92  // if (counter > maxVtx_) break ;
93  // sumpt2 = PVCluster.pTSquaredSum(*vtx) ;
94  // if (sumpt2 > sumpt2previous*fractionSumPt2_ && sumpt2 > minSumPt2_ ) vtxs_trim->push_back(*vtx) ;
95  // else if (counter == 0 ) vtxs_trim->push_back(*vtx) ;
96  // sumpt2previous = sumpt2 ;
97  // }
98 
99  double sumpt2first = pvComparer_->pTSquaredSum(*(vtxs->begin()));
100 
101  for (reco::VertexCollection::const_iterator vtx = vtxs->begin(), evtx = vtxs->end(); vtx != evtx; ++vtx) {
102  if (vtxs_trim->size() >= maxVtx_)
103  break;
104  sumpt2 = pvComparer_->pTSquaredSum(*vtx);
105  // std::cout << "sumpt2: " << sumpt2 << "[" << sumpt2first << "]" << std::endl;
106  // if (sumpt2 >= sumpt2first*fractionSumPt2_ && sumpt2 > minSumPt2_ ) vtxs_trim->push_back(*vtx) ;
107  if (sumpt2 >= sumpt2first * fractionSumPt2_ && sumpt2 > minSumPt2_)
108  vtxs_trim->push_back(*vtx);
109  }
110  // std::cout << " ==> # vertices: " << vtxs_trim->size() << std::endl;
111  iEvent.put(std::move(vtxs_trim));
112 }
113 
114 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
116  //The following says we do not know what parameters are allowed so do no validation
117  // Please change this to state exactly what you do use, even if it is no parameters
119  desc.add<edm::InputTag>("src", edm::InputTag(""))->setComment("input (pixel) vertex collection");
120  desc.add<unsigned int>("maxVtx", 100)->setComment("max output collection size (number of accepted vertices)");
121  desc.add<double>("fractionSumPt2", 0.3)->setComment("threshold on sumPt2 fraction of the leading vertex");
122  desc.add<double>("minSumPt2", 0.)->setComment("min sumPt2");
123  edm::ParameterSetDescription PVcomparerPSet;
124  PVcomparerPSet.add<double>("track_pt_min", 1.0)->setComment("min track p_T");
125  PVcomparerPSet.add<double>("track_pt_max", 10.0)->setComment("max track p_T");
126  PVcomparerPSet.add<double>("track_chi2_max", 99999.)->setComment("max track chi2");
127  PVcomparerPSet.add<double>("track_prob_min", -1.)->setComment("min track prob");
128  desc.add<edm::ParameterSetDescription>("PVcomparer", PVcomparerPSet)
129  ->setComment("from RecoPixelVertexing/PixelVertexFinding/python/PVClusterComparer_cfi.py");
130  descriptions.add("hltPixelVertexCollectionTrimmer", desc);
131 }
132 
133 //define this as a plug-in
ConfigurationDescriptions.h
HLT_FULL_cff.track_pt_min
track_pt_min
Definition: HLT_FULL_cff.py:256
HLT_FULL_cff.track_prob_min
track_prob_min
Definition: HLT_FULL_cff.py:255
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
edm::EDGetTokenT< reco::VertexCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
EDProducer.h
PixelVertexCollectionTrimmer::maxVtx_
unsigned int maxVtx_
Definition: PixelVertexCollectionTrimmer.cc:50
PixelVertexCollectionTrimmer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PixelVertexCollectionTrimmer.cc:78
edm::Handle< reco::VertexCollection >
PixelVertexCollectionTrimmer::PixelVertexCollectionTrimmer
PixelVertexCollectionTrimmer(const edm::ParameterSet &)
Definition: PixelVertexCollectionTrimmer.cc:57
PixelVertexCollectionTrimmer::minSumPt2_
double minSumPt2_
Definition: PixelVertexCollectionTrimmer.cc:52
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
PVClusterComparer::pTSquaredSum
double pTSquaredSum(const PVCluster &v)
Calculate sum of square of the pT's of the tracks in the vertex.
Definition: PVClusterComparer.cc:23
PixelVertexCollectionTrimmer::pvComparer_
PVClusterComparer * pvComparer_
Definition: PixelVertexCollectionTrimmer.cc:54
ParameterSetDescription.h
HLT_FULL_cff.track_chi2_max
track_chi2_max
Definition: HLT_FULL_cff.py:253
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
Vertex.h
HLT_FULL_cff.track_pt_max
track_pt_max
Definition: HLT_FULL_cff.py:254
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:38
PixelVertexCollectionTrimmer::~PixelVertexCollectionTrimmer
~PixelVertexCollectionTrimmer() override
Definition: PixelVertexCollectionTrimmer.cc:75
edm::EventSetup
Definition: EventSetup.h:57
InputTag.h
PixelVertexCollectionTrimmer
Definition: PixelVertexCollectionTrimmer.cc:39
VertexFwd.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
PVClusterComparer.h
PixelVertexCollectionTrimmer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PixelVertexCollectionTrimmer.cc:115
extraflags_cff.vtx
vtx
Definition: extraflags_cff.py:18
Frameworkfwd.h
L1TMuonDQMOffline_cfi.vtxInputTag
vtxInputTag
Definition: L1TMuonDQMOffline_cfi.py:46
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ParameterSet.h
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
PVClusterComparer
Definition: PVClusterComparer.h:16
edm::Event
Definition: Event.h:73
PixelVertexCollectionTrimmer::fractionSumPt2_
double fractionSumPt2_
Definition: PixelVertexCollectionTrimmer.cc:51
edm::InputTag
Definition: InputTag.h:15
PixelVertexCollectionTrimmer::vtxToken_
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
Definition: PixelVertexCollectionTrimmer.cc:49