CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PVSorterBySumP.cc
Go to the documentation of this file.
1 // Includes
8 
10 
13 
14 #include <memory>
15 #include <vector>
16 #include <sstream>
17 
18 
20 // class definition
22 template<typename T1>
24 {
25 public:
26  // construction/destruction
27  bestPVselector(const edm::ParameterSet& iConfig);
28  virtual ~bestPVselector();
29 
30  // member functions
31  void produce(edm::Event& iEvent,const edm::EventSetup& iSetup) override;
32  void endJob() override;
33 
34 private:
35  // member data
37 };
38 
39 
40 
42 // construction/destruction
44 
45 //______________________________________________________________________________
46 template<typename T1>
48  : src_(consumes<std::vector<T1> >(iConfig.getParameter<edm::InputTag>("src")))
49 {
50  produces<std::vector<T1> >();
51 }
52 
53 
54 //______________________________________________________________________________
55 template<typename T1>
57 
59 // implementation of member functions
61 
62 //______________________________________________________________________________
63 template<typename T1>
65 {
66  std::auto_ptr<std::vector<T1> > theBestPV(new std::vector<T1 >);
67 
68  edm::Handle< std::vector<T1> > VertexHandle;
69  iEvent.getByToken(src_,VertexHandle);
70 
71  if( VertexHandle->size() == 0 )
72  {
73  iEvent.put(theBestPV);
74  return ;
75  }
76 
77  typename std::vector<T1>::const_iterator PVit ;
78  typename std::vector<T1>::const_iterator bestPV ;
79 
80  double bestP4 = 0 ;
81  double sumSquarePt = 0 ;
82 
83  for (PVit = VertexHandle->begin(); PVit != VertexHandle->end(); ++PVit) {
84  sumSquarePt = (PVit -> p4().pt())*(PVit -> p4().pt()) ;
85  if( sumSquarePt > bestP4 ){
86  bestP4 = sumSquarePt ;
87  bestPV = PVit ;
88  }
89  }
90 
91  theBestPV->push_back( *bestPV );
92  iEvent.put(theBestPV);
93 
94 }
95 
96 template<typename T1>
98 {
99 }
100 
102 
bestPVselector< reco::Vertex > HighestSumP4PrimaryVertexSelector
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bestPVselector(const edm::ParameterSet &iConfig)
return((rh^lh)&mask)
void endJob() override
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
double p4[4]
Definition: TauolaWrapper.h:92
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
virtual ~bestPVselector()
edm::EDGetTokenT< std::vector< T1 > > src_