CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EventVtxInfoNtupleDumper.cc
Go to the documentation of this file.
8 
9 #include <vector>
10 
11 using namespace edm;
12 using namespace std;
13 using namespace reco;
14 
15 
17 public:
19 
20 private:
21  void produce( edm::Event &, const edm::EventSetup & );
23 
24 };
25 
27  primaryVertices_(cfg.getParameter<InputTag>("primaryVertices")) {
28  produces<int>( "numPV" ).setBranchAlias( "numPV" );
29  produces<int>( "nTrkPV" ).setBranchAlias( "nTrkPV" );
30  produces<float>( "chi2PV" ).setBranchAlias( "chi2PV" );
31  produces<float>( "ndofPV" ).setBranchAlias( "ndofPV" );
32  produces<float>( "zPV" ).setBranchAlias( "zPV" );
33  produces<float>( "rhoPV" ).setBranchAlias( "rhoPV" );
34  // produces<std::vector< unsigned int > >( "nTrkPV" ).setBranchAlias( "nTrkPV" );
35  // produces<std::vector< float > >( "chi2PV" ).setBranchAlias( "chi2PV" );
36  // produces<std::vector< float > >( "ndofPV" ).setBranchAlias( "ndofPV" );
37 }
38 
39 
40 
42 
43  Handle<reco::VertexCollection> primaryVertices; // Collection of primary Vertices
44  evt.getByLabel(primaryVertices_, primaryVertices);
45  auto_ptr<int> nVtxs( new int );
46  auto_ptr<int> nTrkVtx( new int );
47  auto_ptr<float> chi2Vtx( new float );
48  auto_ptr<float> ndofVtx( new float );
49  auto_ptr<float> zVtx( new float );
50  auto_ptr<float> rhoVtx( new float );
51  // auto_ptr< vector< unsigned int > > nTrkVtx( new vector< unsigned int > );
52  // auto_ptr< vector< float > > chi2Vtx( new vector< float > );
53  // auto_ptr< vector< float > > ndofVtx( new vector< float > );
54 
55  const reco::Vertex &pv = (*primaryVertices)[0];
56 
57  *nVtxs = -1;
58  *nTrkVtx = -1;
59  *chi2Vtx = -1.0;
60  *ndofVtx = -1.0;
61  *zVtx = -1000;
62  *rhoVtx = -1000;
63  if( !(pv.isFake()) ) {
64  *nVtxs = primaryVertices->size();
65  *nTrkVtx = pv.tracksSize();
66  *chi2Vtx = pv.chi2();
67  *ndofVtx = pv.ndof();
68  *zVtx = pv.z();
69  *rhoVtx = pv.position().Rho();
70  }
71  // nTrkVtx->push_back(pv.tracksSize());
72  // chi2Vtx->push_back(pv.chi2());
73  // ndofVtx->push_back(pv.ndof());
74  evt.put( nVtxs, "numPV" );
75  evt.put( nTrkVtx, "nTrkPV" );
76  evt.put( chi2Vtx, "chi2PV" );
77  evt.put( ndofVtx, "ndofPV" );
78  evt.put( zVtx, "zPV" );
79  evt.put( rhoVtx, "rhoPV" );
80 
81 }
82 
84 
86 
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const Point & position() const
position
Definition: Vertex.h:93
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
double chi2() const
chi-squares
Definition: Vertex.h:82
double z() const
y coordinate
Definition: Vertex.h:99
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double ndof() const
Definition: Vertex.h:89
EventVtxInfoNtupleDumper(const edm::ParameterSet &)
bool isFake() const
Definition: Vertex.h:65
void produce(edm::Event &, const edm::EventSetup &)
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:35