CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoHI/HiTracking/plugins/HIPixelClusterVtxProducer.cc

Go to the documentation of this file.
00001 #include "HIPixelClusterVtxProducer.h"
00002 
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "FWCore/Framework/interface/Frameworkfwd.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "DataFormats/VertexReco/interface/Vertex.h"
00009 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00010 
00011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00012 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00013 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00014 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00015 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00016 
00017 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00018 
00019 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00020 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00021 
00022 
00023 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00024 
00025 #include <fstream>
00026 #include <iostream>
00027 #include <vector>
00028 #include <algorithm>
00029 
00030 /*****************************************************************************/
00031 HIPixelClusterVtxProducer::HIPixelClusterVtxProducer(const edm::ParameterSet& ps)
00032   : srcPixels_(ps.getParameter<std::string>("pixelRecHits")),
00033     minZ_(ps.getParameter<double>("minZ")),
00034     maxZ_(ps.getParameter<double>("maxZ")),
00035     zStep_(ps.getParameter<double>("zStep"))
00036 
00037 {
00038   // Constructor
00039   produces<reco::VertexCollection>();
00040 }
00041 
00042 
00043 /*****************************************************************************/
00044 HIPixelClusterVtxProducer::~HIPixelClusterVtxProducer()
00045 { 
00046   // Destructor
00047 }
00048 
00049 
00050 /*****************************************************************************/
00051 void HIPixelClusterVtxProducer::produce(edm::Event& ev, const edm::EventSetup& es)
00052 {
00053 
00054   // new vertex collection
00055   std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
00056 
00057   // get pixel rechits
00058   edm::Handle<SiPixelRecHitCollection> hRecHits;
00059   ev.getByLabel(edm::InputTag(srcPixels_),hRecHits);
00060 
00061   // get tracker geometry
00062   if (hRecHits.isValid()) {
00063     edm::ESHandle<TrackerGeometry> trackerHandle;
00064     es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
00065     const TrackerGeometry *tgeo = trackerHandle.product();
00066     const SiPixelRecHitCollection *hits = hRecHits.product();
00067 
00068     // loop over pixel rechits
00069     std::vector<VertexHit> vhits;
00070     for(SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(), 
00071           end = hits->data().end(); hit != end; ++hit) {
00072       if (!hit->isValid())
00073         continue;
00074       DetId id(hit->geographicalId());
00075       if(id.subdetId() != int(PixelSubdetector::PixelBarrel))
00076         continue;
00077       const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit*>(tgeo->idToDet(id));
00078       if (1) {
00079         const PixelTopology *pixTopo = &(pgdu->specificTopology());
00080         std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
00081         bool pixelOnEdge = false;
00082         for(std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); 
00083             pixel != pixels.end(); ++pixel) {
00084           int pixelX = pixel->x;
00085           int pixelY = pixel->y;
00086           if(pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
00087             pixelOnEdge = true;
00088             break;
00089           }
00090         }
00091         if (pixelOnEdge)
00092           continue;
00093       }
00094 
00095       LocalPoint lpos = LocalPoint(hit->localPosition().x(),
00096                                    hit->localPosition().y(),
00097                                    hit->localPosition().z());
00098       GlobalPoint gpos = pgdu->toGlobal(lpos);
00099       VertexHit vh;
00100       vh.z = gpos.z(); 
00101       vh.r = gpos.perp(); 
00102       vh.w = hit->cluster()->sizeY();
00103       vhits.push_back(vh);
00104     }
00105     
00106     // estimate z-position from cluster lengths
00107     double zest = 0.0;
00108     int nhits = 0, nhits_max = 0;
00109     double chi = 0, chi_max = 1e+9;
00110     for(double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
00111       nhits = getContainedHits(vhits, z0, chi);
00112       if(nhits == 0) 
00113         continue;
00114       if(nhits > nhits_max) { 
00115         chi_max = 1e+9; 
00116         nhits_max = nhits; 
00117       }
00118       if(nhits >= nhits_max && chi < chi_max) {
00119         chi_max = chi; 
00120         zest = z0;   
00121       }
00122     } 
00123 
00124     LogTrace("MinBiasTracking")
00125       << "  [vertex position] estimated = " << zest 
00126       << " | pixel barrel hits = " << vhits.size();
00127 
00128     // put 1-d vertex and dummy errors into collection
00129     reco::Vertex::Error err;
00130     err(2,2) = 0.6 * 0.6;
00131     reco::Vertex ver(reco::Vertex::Point(0,0,zest), err, 0, 1, 1);
00132     vertices->push_back(ver);
00133   }
00134 
00135   ev.put(vertices);
00136 }
00137 
00138 /*****************************************************************************/
00139 int HIPixelClusterVtxProducer::getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi)
00140 {
00141   // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
00142   int n = 0;
00143   chi   = 0.;
00144 
00145   for(std::vector<VertexHit>::const_iterator hit = hits.begin(); hit!= hits.end(); hit++) {
00146     double p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
00147     if(TMath::Abs(p - hit->w) <= 1.) { 
00148       chi += fabs(p - hit->w);
00149       n++;
00150     }
00151   }
00152   return n;
00153 }