Go to the documentation of this file.00001 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
00002
00003
00004
00005
00006
00007 class HLTPixelClusterShapeFilter : public HLTFilter {
00008 public:
00009 explicit HLTPixelClusterShapeFilter(const edm::ParameterSet&);
00010 ~HLTPixelClusterShapeFilter();
00011
00012 private:
00013
00014 edm::InputTag inputTag_;
00015 bool saveTag_;
00016
00017 double minZ_;
00018 double maxZ_;
00019 double zStep_;
00020
00021 std::vector<double> clusterPars_;
00022 int nhitsTrunc_;
00023 double clusterTrunc_;
00024
00025 struct VertexHit
00026 {
00027 float z;
00028 float r;
00029 float w;
00030 };
00031
00032 virtual bool filter(edm::Event&, const edm::EventSetup&);
00033 int getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi);
00034
00035 };
00036
00037 #include "FWCore/Framework/interface/Event.h"
00038 #include "FWCore/Framework/interface/EventSetup.h"
00039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00040 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00041 #include "DataFormats/Common/interface/Handle.h"
00042 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00043 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00044
00045 #include "FWCore/Framework/interface/ESHandle.h"
00046 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00047 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00048 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00049 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00050 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00051 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00052 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00053 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00054 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00055
00056
00057
00058
00059
00060 HLTPixelClusterShapeFilter::HLTPixelClusterShapeFilter(const edm::ParameterSet& config) :
00061 inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
00062 saveTag_ (config.getUntrackedParameter<bool>("saveTag", false)),
00063 minZ_ (config.getParameter<double>("minZ")),
00064 maxZ_ (config.getParameter<double>("maxZ")),
00065 zStep_ (config.getParameter<double>("zStep")),
00066 clusterPars_ (config.getParameter< std::vector<double> >("clusterPars")),
00067 nhitsTrunc_ (config.getParameter<int>("nhitsTrunc")),
00068 clusterTrunc_ (config.getParameter<double>("clusterTrunc"))
00069 {
00070 LogDebug("") << "Using the " << inputTag_ << " input collection";
00071
00072
00073 produces<trigger::TriggerFilterObjectWithRefs>();
00074 }
00075
00076 HLTPixelClusterShapeFilter::~HLTPixelClusterShapeFilter()
00077 {
00078 }
00079
00080
00081
00082
00083
00084
00085 bool HLTPixelClusterShapeFilter::filter(edm::Event& event, const edm::EventSetup& iSetup)
00086 {
00087
00088
00089
00090
00091
00092 std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterobject (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00093 if (saveTag_) filterobject->addCollectionTag(inputTag_);
00094 bool accept = true;
00095
00096
00097 edm::Handle<SiPixelRecHitCollection> hRecHits;
00098 event.getByLabel(inputTag_, hRecHits);
00099
00100
00101 if (hRecHits.isValid()) {
00102 edm::ESHandle<TrackerGeometry> trackerHandle;
00103 iSetup.get<TrackerDigiGeometryRecord>().get(trackerHandle);
00104 const TrackerGeometry *tgeo = trackerHandle.product();
00105 const SiPixelRecHitCollection *hits = hRecHits.product();
00106
00107
00108 int nPxlHits=0;
00109 std::vector<VertexHit> vhits;
00110 for(SiPixelRecHitCollection::DataContainer::const_iterator hit = hits->data().begin(),
00111 end = hits->data().end(); hit != end; ++hit) {
00112 if (!hit->isValid())
00113 continue;
00114 ++nPxlHits;
00115 DetId id(hit->geographicalId());
00116 if(id.subdetId() != int(PixelSubdetector::PixelBarrel))
00117 continue;
00118 const PixelGeomDetUnit *pgdu = static_cast<const PixelGeomDetUnit*>(tgeo->idToDet(id));
00119 if (1) {
00120 const PixelTopology *pixTopo = &(pgdu->specificTopology());
00121 std::vector<SiPixelCluster::Pixel> pixels(hit->cluster()->pixels());
00122 bool pixelOnEdge = false;
00123 for(std::vector<SiPixelCluster::Pixel>::const_iterator pixel = pixels.begin();
00124 pixel != pixels.end(); ++pixel) {
00125 int pixelX = pixel->x;
00126 int pixelY = pixel->y;
00127 if(pixTopo->isItEdgePixelInX(pixelX) || pixTopo->isItEdgePixelInY(pixelY)) {
00128 pixelOnEdge = true;
00129 break;
00130 }
00131 }
00132 if (pixelOnEdge)
00133 continue;
00134 }
00135
00136 LocalPoint lpos = LocalPoint(hit->localPosition().x(),
00137 hit->localPosition().y(),
00138 hit->localPosition().z());
00139 GlobalPoint gpos = pgdu->toGlobal(lpos);
00140 VertexHit vh;
00141 vh.z = gpos.z();
00142 vh.r = gpos.perp();
00143 vh.w = hit->cluster()->sizeY();
00144 vhits.push_back(vh);
00145 }
00146
00147
00148 double zest = 0.0;
00149 int nhits = 0, nhits_max = 0;
00150 double chi = 0, chi_max = 1e+9;
00151 for(double z0 = minZ_; z0 <= maxZ_; z0 += zStep_) {
00152 nhits = getContainedHits(vhits, z0, chi);
00153 if(nhits == 0)
00154 continue;
00155 if(nhits > nhits_max) {
00156 chi_max = 1e+9;
00157 nhits_max = nhits;
00158 }
00159 if(nhits >= nhits_max && chi < chi_max) {
00160 chi_max = chi;
00161 zest = z0;
00162 }
00163 }
00164
00165 chi = 0;
00166 int nbest=0, nminus=0, nplus=0;
00167 nbest = getContainedHits(vhits,zest,chi);
00168 nminus = getContainedHits(vhits,zest-10.,chi);
00169 nplus = getContainedHits(vhits,zest+10.,chi);
00170
00171 double clusVtxQual=0.0;
00172 if ((nminus+nplus)> 0)
00173 clusVtxQual = (2.0*nbest)/(nminus+nplus);
00174 else if (nbest>0)
00175 clusVtxQual = 1000.0;
00176 else
00177 clusVtxQual = 0;
00178
00179
00180 double polyCut=0;
00181 for(unsigned int i=0; i < clusterPars_.size(); i++) {
00182 polyCut += clusterPars_[i]*std::pow((double)nPxlHits,(int)i);
00183 }
00184 if(nPxlHits < nhitsTrunc_)
00185 polyCut=0;
00186 if(polyCut > clusterTrunc_ && clusterTrunc_ > 0)
00187 polyCut=clusterTrunc_;
00188
00189 if (clusVtxQual < polyCut) accept = false;
00190
00191
00192 event.put(filterobject);
00193 }
00194
00195
00196 return accept;
00197 }
00198
00199 int HLTPixelClusterShapeFilter::getContainedHits(const std::vector<VertexHit> &hits, double z0, double &chi)
00200 {
00201
00202 int n = 0;
00203 chi = 0.;
00204
00205 for(std::vector<VertexHit>::const_iterator hit = hits.begin(); hit!= hits.end(); hit++) {
00206 double p = 2 * fabs(hit->z - z0)/hit->r + 0.5;
00207 if(fabs(p - hit->w) <= 1.) {
00208 chi += fabs(p - hit->w);
00209 n++;
00210 }
00211 }
00212 return n;
00213 }
00214
00215
00216
00217 #include "FWCore/Framework/interface/MakerMacros.h"
00218 DEFINE_FWK_MODULE(HLTPixelClusterShapeFilter);