CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

HIPixelClusterVtxProducer Class Reference

#include <HIPixelClusterVtxProducer.h>

Inheritance diagram for HIPixelClusterVtxProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Classes

struct  VertexHit

Public Member Functions

 HIPixelClusterVtxProducer (const edm::ParameterSet &ps)
 ~HIPixelClusterVtxProducer ()

Private Member Functions

int getContainedHits (const std::vector< VertexHit > &hits, double z0, double &chi)
virtual void produce (edm::Event &ev, const edm::EventSetup &es)

Private Attributes

double maxZ_
double minZ_
std::string srcPixels_
double zStep_

Detailed Description

Definition at line 11 of file HIPixelClusterVtxProducer.h.


Constructor & Destructor Documentation

HIPixelClusterVtxProducer::HIPixelClusterVtxProducer ( const edm::ParameterSet ps) [explicit]

Definition at line 32 of file HIPixelClusterVtxProducer.cc.

  : srcPixels_(ps.getParameter<std::string>("pixelRecHits")),
    minZ_(ps.getParameter<double>("minZ")),
    maxZ_(ps.getParameter<double>("maxZ")),
    zStep_(ps.getParameter<double>("zStep"))

{
  // Constructor
  produces<reco::VertexCollection>();
}
HIPixelClusterVtxProducer::~HIPixelClusterVtxProducer ( )

Definition at line 45 of file HIPixelClusterVtxProducer.cc.

{ 
  // Destructor
}

Member Function Documentation

int HIPixelClusterVtxProducer::getContainedHits ( const std::vector< VertexHit > &  hits,
double  z0,
double &  chi 
) [private]

Definition at line 142 of file HIPixelClusterVtxProducer.cc.

References n, AlCaHLTBitMon_ParallelJobs::p, and hit::z.

Referenced by produce().

{
  // Calculate number of hits contained in v-shaped window in cluster y-width vs. z-position.
  int n = 0;
  chi   = 0.;

  for(std::vector<VertexHit>::const_iterator hit = hits.begin(); hit!= hits.end(); hit++) {
    double p = 2 * fabs(hit->z - z0)/hit->r + 0.5; // FIXME
    if(TMath::Abs(p - hit->w) <= 1.) { 
      chi += fabs(p - hit->w);
      n++;
    }
  }
  return n;
}
void HIPixelClusterVtxProducer::produce ( edm::Event ev,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 52 of file HIPixelClusterVtxProducer.cc.

References edmNew::DetSetVector< T >::data(), alignCSCRings::e, end, edm::EventSetup::get(), edm::Event::getByLabel(), getContainedHits(), TrackerGeometry::idToDet(), PixelTopology::isItEdgePixelInX(), PixelTopology::isItEdgePixelInY(), edm::HandleBase::isValid(), LogTrace, maxZ_, minZ_, PixelSubdetector::PixelBarrel, edm::ESHandle< T >::product(), edm::Handle< T >::product(), edm::Event::put(), HIPixelClusterVtxProducer::VertexHit::r, PixelGeomDetUnit::specificTopology(), srcPixels_, GeomDet::toGlobal(), cmsDownloadME::ver, HIPixelClusterVtxProducer::VertexHit::w, hit::x, hit::y, HIPixelClusterVtxProducer::VertexHit::z, hit::z, and zStep_.

{

  // new vertex collection
  std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);

  // get pixel rechits
  edm::Handle<SiPixelRecHitCollection> hRecHits;
  try {
    ev.getByLabel(edm::InputTag(srcPixels_),hRecHits);
  } catch (...) {}

  // get tracker geometry
  if (hRecHits.isValid()) {
    edm::ESHandle<TrackerGeometry> trackerHandle;
    es.get<TrackerDigiGeometryRecord>().get(trackerHandle);
    const TrackerGeometry *tgeo = trackerHandle.product();
    const SiPixelRecHitCollection *hits = hRecHits.product();

    // loop over pixel rechits
    std::vector<VertexHit> vhits;
    for(SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(), 
          end = hits->data().end(); hit != end; ++hit) {
      if (!hit->isValid())
        continue;
      DetId id(hit->geographicalId());
      if(id.subdetId() != int(PixelSubdetector::PixelBarrel))
        continue;
      const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit*>(tgeo->idToDet(id));
      if (1) {
        const PixelTopology *pixTopo = &(pgdu->specificTopology());
        std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
        bool pixelOnEdge = false;
        for(std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin(); 
            pixel != pixels.end(); ++pixel) {
          int pixelX = pixel->x;
          int pixelY = pixel->y;
          if(pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
            pixelOnEdge = true;
            break;
          }
        }
        if (pixelOnEdge)
          continue;
      }

      LocalPoint lpos = LocalPoint(hit->localPosition().x(),
                                   hit->localPosition().y(),
                                   hit->localPosition().z());
      GlobalPoint gpos = pgdu->toGlobal(lpos);
      VertexHit vh;
      vh.z = gpos.z(); 
      vh.r = gpos.perp(); 
      vh.w = hit->cluster()->sizeY();
      vhits.push_back(vh);
    }
    
    // estimate z-position from cluster lengths
    double zest = 0.0;
    int nhits = 0, nhits_max = 0;
    double chi = 0, chi_max = 1e+9;
    for(double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
      nhits = getContainedHits(vhits, z0, chi);
      if(nhits == 0) 
        continue;
      if(nhits > nhits_max) { 
        chi_max = 1e+9; 
        nhits_max = nhits; 
      }
      if(nhits >= nhits_max && chi < chi_max) {
        chi_max = chi; 
        zest = z0;   
      }
    } 

    LogTrace("MinBiasTracking")
      << "  [vertex position] estimated = " << zest 
      << " | pixel barrel hits = " << vhits.size();

    // put 1-d vertex and dummy errors into collection
    reco::Vertex::Error err;
    err(2,2) = 0.6 * 0.6;
    reco::Vertex ver(reco::Vertex::Point(0,0,zest), err, 0, 1, 1);
    vertices->push_back(ver);
  }

  ev.put(vertices);
}

Member Data Documentation

Definition at line 31 of file HIPixelClusterVtxProducer.h.

Referenced by produce().

Definition at line 30 of file HIPixelClusterVtxProducer.h.

Referenced by produce().

Definition at line 28 of file HIPixelClusterVtxProducer.h.

Referenced by produce().

Definition at line 32 of file HIPixelClusterVtxProducer.h.

Referenced by produce().