CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PixelVertexProducer Class Reference

#include <PixelVertexFinding/interface/PixelVertexProducer.h>

Inheritance diagram for PixelVertexProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 PixelVertexProducer (const edm::ParameterSet &)
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
 ~PixelVertexProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

DivisiveVertexFinderdvf_
 
bool method2
 
double ptMin_
 
edm::EDGetTokenT< reco::BeamSpottoken_BeamSpot
 
edm::EDGetTokenT
< reco::TrackCollection
token_Tracks
 
edm::InputTag trackCollName
 
int verbose_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: This produces 1D (z only) primary vertexes using only pixel information.

Implementation: This producer can use either the Divisive Primary Vertex Finder or the Histogramming Primary Vertex Finder (currently not implemented). It relies on the PixelTripletProducer and PixelTrackProducer having already been run upstream. This is code ported from ORCA originally written by S Cucciarelli, M Konecki, D Kotlinski.

Definition at line 35 of file PixelVertexProducer.h.

Constructor & Destructor Documentation

PixelVertexProducer::PixelVertexProducer ( const edm::ParameterSet conf)
explicit

Definition at line 13 of file PixelVertexProducer.cc.

References dvf_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), method2, ptMin_, AlCaHLTBitMon_QueryRunRegistry::string, token_BeamSpot, token_Tracks, trackCollName, and verbose_.

14  : verbose_(0), dvf_(0), ptMin_(1.0)
15 {
16  // Register my product
17  produces<reco::VertexCollection>();
18 
19  // Setup shop
20  verbose_ = conf.getParameter<int>("Verbosity"); // 0 silent, 1 chatty, 2 loud
21  std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder
22  bool useError = conf.getParameter<bool>("UseError"); // true
23  bool wtAverage = conf.getParameter<bool>("WtAverage"); // true
24  double zOffset = conf.getParameter<double>("ZOffset"); // 5.0 sigma
25  double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm
26  int ntrkMin = conf.getParameter<int>("NTrkMin"); // 3
27  // Tracking requirements before sending a track to be considered for vtx
28  ptMin_ = conf.getParameter<double>("PtMin"); // 1.0 GeV
29  trackCollName = conf.getParameter<edm::InputTag>("TrackCollection");
30  token_Tracks = consumes<reco::TrackCollection>(trackCollName);
31  token_BeamSpot = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpot"));
32  method2 = conf.getParameter<bool>("Method2");
33 
34  double track_pt_min = ptMin_;
35  double track_pt_max = 10.;
36  double track_chi2_max = 9999999.;
37  double track_prob_min = -1.;
38 
39  if ( conf.exists("PVcomparer") ) {
40  edm::ParameterSet PVcomparerPSet = conf.getParameter<edm::ParameterSet>("PVcomparer");
41  track_pt_min = PVcomparerPSet.getParameter<double>("track_pt_min");
42  if (track_pt_min != ptMin_) {
43  if (track_pt_min < ptMin_)
44  edm::LogWarning("PixelVertexProducer") << "minimum track pT setting differs between PixelVertexProducer (" << ptMin_ << ") and PVcomparer (" << track_pt_min << ") [PVcomparer considers tracks w/ lower threshold than PixelVertexProducer does] !!!";
45  else
46  edm::LogInfo("PixelVertexProducer") << "minimum track pT setting differs between PixelVertexProducer (" << ptMin_ << ") and PVcomparer (" << track_pt_min << ") !!!";
47  }
48  track_pt_max = PVcomparerPSet.getParameter<double>("track_pt_max");
49  track_chi2_max = PVcomparerPSet.getParameter<double>("track_chi2_max");
50  track_prob_min = PVcomparerPSet.getParameter<double>("track_prob_min");
51  }
52 
53  if (finder == "DivisiveVertexFinder") {
54  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n";
55  dvf_ = new DivisiveVertexFinder(track_pt_min,track_pt_max,track_chi2_max,track_prob_min,zOffset, ntrkMin, useError, zSeparation, wtAverage, verbose_);
56  }
57  else { // Finder not supported, or you made a mistake in your request
58  // throw an exception once I figure out how CMSSW does this
59  }
60 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::InputTag trackCollName
edm::EDGetTokenT< reco::TrackCollection > token_Tracks
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot
DivisiveVertexFinder * dvf_
PixelVertexProducer::~PixelVertexProducer ( )

Definition at line 63 of file PixelVertexProducer.cc.

References dvf_.

63  {
64  delete dvf_;
65 }
DivisiveVertexFinder * dvf_

Member Function Documentation

void PixelVertexProducer::produce ( edm::Event e,
const edm::EventSetup es 
)
virtual

Implements edm::EDProducer.

Definition at line 67 of file PixelVertexProducer.cc.

References reco::Vertex::add(), dvf_, reco::BeamSpot::dxdz(), reco::BeamSpot::dydz(), relativeConstraints::error, DivisiveVertexFinder::findVertexes(), DivisiveVertexFinder::findVertexesAlt(), edm::Event::getByToken(), i, edm::HandleBase::isValid(), method2, convertSQLiteXML::ok, reco::BeamSpot::position(), edm::Handle< T >::product(), EnergyCorrector::pt, ptMin_, edm::RefVector< C, T, F >::push_back(), edm::Event::put(), reco::BeamSpot::rotatedCovariance3D(), edm::RefVector< C, T, F >::size(), mathSSE::sqrt(), token_BeamSpot, token_Tracks, trackCollName, testEve_cfg::tracks, findQualityFiles::v, verbose_, x, reco::BeamSpot::x0(), detailsBasic3DVector::y, reco::BeamSpot::y0(), detailsBasic3DVector::z, and reco::BeamSpot::z0().

67  {
68 
69  // First fish the pixel tracks out of the event
70  edm::Handle<reco::TrackCollection> trackCollection;
71  e.getByToken(token_Tracks,trackCollection);
72  const reco::TrackCollection tracks = *(trackCollection.product());
73  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called " << trackCollName << "\n";
74 
75 
76  // Second, make a collection of pointers to the tracks we want for the vertex finder
78  for (unsigned int i=0; i<tracks.size(); i++) {
79  if (tracks[i].pt() > ptMin_)
80  trks.push_back( reco::TrackRef(trackCollection, i) );
81  }
82  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n";
83 
85  e.getByToken(token_BeamSpot,bsHandle);
86  math::XYZPoint myPoint(0.,0.,0.);
87  if (bsHandle.isValid()) myPoint = math::XYZPoint(bsHandle->x0(),bsHandle->y0(), 0. ); //FIXME: fix last coordinate with vertex.z() at same time
88 
89  // Third, ship these tracks off to be vertexed
90  std::auto_ptr<reco::VertexCollection> vertexes(new reco::VertexCollection);
91  bool ok;
92  if (method2) {
93  ok = dvf_->findVertexesAlt(trks, // input
94  *vertexes,myPoint); // output
95  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok;
96  }
97  else {
98  ok = dvf_->findVertexes(trks, // input
99  *vertexes); // output
100  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok;
101  }
102 
103  if (verbose_ > 0) {
104  edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n";
105  for (unsigned int i=0; i<vertexes->size(); ++i) {
106  edm::LogInfo("PixelVertexProducer") << "Vertex number " << i << " has " << (*vertexes)[i].tracksSize() << " tracks with a position of " << (*vertexes)[i].z() << " +- " << std::sqrt( (*vertexes)[i].covariance(2,2) );
107  }
108  }
109 
110 
111  if(bsHandle.isValid())
112  {
113  const reco::BeamSpot & bs = *bsHandle;
114 
115  for (unsigned int i=0; i<vertexes->size(); ++i) {
116  double z=(*vertexes)[i].z();
117  double x=bs.x0()+bs.dxdz()*(z-bs.z0());
118  double y=bs.y0()+bs.dydz()*(z-bs.z0());
119  reco::Vertex v( reco::Vertex::Point(x,y,z), (*vertexes)[i].error(),(*vertexes)[i].chi2() , (*vertexes)[i].ndof() , (*vertexes)[i].tracksSize());
120  //Copy also the tracks
121  for (std::vector<reco::TrackBaseRef >::const_iterator it = (*vertexes)[i].tracks_begin();
122  it !=(*vertexes)[i].tracks_end(); it++) {
123  v.add( *it );
124  }
125  (*vertexes)[i]=v;
126 
127  }
128  }
129  else
130  {
131  edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) ";
132  }
133 
134  if(!vertexes->size() && bsHandle.isValid()){
135 
136  const reco::BeamSpot & bs = *bsHandle;
137 
139  if ( (bse.cxx() <= 0.) ||
140  (bse.cyy() <= 0.) ||
141  (bse.czz() <= 0.) ) {
143  we(0,0)=10000;
144  we(1,1)=10000;
145  we(2,2)=10000;
146  vertexes->push_back(reco::Vertex(bs.position(), we,0.,0.,0));
147 
148  edm::LogInfo("PixelVertexProducer") <<"No vertices found. Beamspot with invalid errors " << bse.matrix() << std::endl
149  << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n"
150  << (*vertexes)[0].x() << "\n"
151  << (*vertexes)[0].y() << "\n"
152  << (*vertexes)[0].z() << "\n";
153  } else {
154  vertexes->push_back(reco::Vertex(bs.position(),
155  bs.rotatedCovariance3D(),0.,0.,0));
156 
157  edm::LogInfo("PixelVertexProducer") << "No vertices found. Will put Vertex derived from BeamSpot into Event:\n"
158  << (*vertexes)[0].x() << "\n"
159  << (*vertexes)[0].y() << "\n"
160  << (*vertexes)[0].z() << "\n";
161  }
162  }
163 
164  else if(!vertexes->size() && !bsHandle.isValid())
165  {
166  edm::LogWarning("PixelVertexProducer") << "No beamspot and no vertex found. No vertex returned.";
167  }
168 
169  e.put(vertexes);
170 
171 }
double z0() const
z coordinate
Definition: BeamSpot.h:68
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
bool findVertexesAlt(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes, const math::XYZPoint &bs)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::InputTag trackCollName
float float float z
edm::EDGetTokenT< reco::TrackCollection > token_Tracks
double dydz() const
dydz slope
Definition: BeamSpot.h:84
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot
T sqrt(T t)
Definition: SSEVec.h:48
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
bool findVertexes(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes)
Run the divisive algorithm and return a vector of vertexes for the input track collection.
bool isValid() const
Definition: HandleBase.h:76
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
DivisiveVertexFinder * dvf_
double y0() const
y coordinate
Definition: BeamSpot.h:66
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
const Point & position() const
position
Definition: BeamSpot.h:62
Definition: DDAxes.h:10
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:78
double x0() const
x coordinate
Definition: BeamSpot.h:64

Member Data Documentation

DivisiveVertexFinder* PixelVertexProducer::dvf_
private

Definition at line 45 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), produce(), and ~PixelVertexProducer().

bool PixelVertexProducer::method2
private

Definition at line 51 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().

double PixelVertexProducer::ptMin_
private

Definition at line 47 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().

edm::EDGetTokenT<reco::BeamSpot> PixelVertexProducer::token_BeamSpot
private

Definition at line 50 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().

edm::EDGetTokenT<reco::TrackCollection> PixelVertexProducer::token_Tracks
private

Definition at line 49 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().

edm::InputTag PixelVertexProducer::trackCollName
private

Definition at line 48 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().

int PixelVertexProducer::verbose_
private

Definition at line 44 of file PixelVertexProducer.h.

Referenced by PixelVertexProducer(), and produce().