CMS 3D CMS Logo

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