CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/EventFilter/SiPixelRawToDigi/src/PixelUnpackingRegions.cc

Go to the documentation of this file.
00001 //
00002 // $Id: PixelUnpackingRegions.cc,v 1.3 2012/03/19 06:58:18 khotilov Exp $
00003 //
00004 #include "EventFilter/SiPixelRawToDigi/interface/PixelUnpackingRegions.h"
00005 
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/ESTransientHandle.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00012 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00013 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00014 
00015 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00016 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00017 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00018 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00019 #include "DataFormats/Math/interface/normalizedPhi.h"
00020 
00021 #include <algorithm>
00022 #include <iterator>
00023 
00024 // local convenience functions
00025 namespace 
00026 {
00027 bool isBPIXFED(unsigned int fed) {return fed< 32;}
00028 bool isFPIXFED(unsigned int fed) {return fed>=32;}
00029 bool isBPIXModule(unsigned int id) {return DetId(id).subdetId() == PixelSubdetector::PixelBarrel;}
00030 bool isFPIXModule(unsigned int id) {return DetId(id).subdetId() == PixelSubdetector::PixelEndcap;}
00031 
00032 std::ostream& operator<<(std::ostream& s, const PixelUnpackingRegions::Module& m)
00033 {
00034   s<< (isBPIXModule(m.id) ? "BPIX " : "FPIX ") <<m.id<<" "<<m.fed<<"   "<<m.phi<<"   "<<m.x<<" "<<m.y<<" "<<m.z<<"  "<<sqrt(std::pow(m.x,2)+std::pow(m.y,2));
00035   return s;
00036 }
00037 }
00038 
00039 
00040 
00041 PixelUnpackingRegions::PixelUnpackingRegions(const edm::ParameterSet& conf)
00042 {
00043   edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("Regions");
00044   beamSpotTag_ = regPSet.getParameter<edm::InputTag>("beamSpot");
00045   inputs_      = regPSet.getParameter<std::vector<edm::InputTag> >("inputs");
00046   dPhi_ = regPSet.getParameter<std::vector<double> >("deltaPhi");
00047   maxZ_ = regPSet.getParameter<std::vector<double> >("maxZ");
00048 
00049   if (inputs_.size() != dPhi_.size() || dPhi_.size() != maxZ_.size() )
00050   {
00051     edm::LogError("PixelUnpackingRegions")<<"Not the same size of config parameters vectors!\n"
00052         <<"   inputs "<<inputs_.size()<<"  deltaPhi "<<dPhi_.size() <<"  maxZ "<< maxZ_.size();
00053   }
00054 }
00055 
00056 
00057 void PixelUnpackingRegions::run(const edm::Event& e, const edm::EventSetup& es)
00058 {
00059   feds_.clear();
00060   modules_.clear();
00061   nreg_ = 0;
00062 
00063   initialize(es);
00064 
00065   edm::Handle<reco::BeamSpot> beamSpot;
00066   e.getByLabel(beamSpotTag_, beamSpot);
00067   beamSpot_ = beamSpot->position();
00068   //beamSpot_ = math::XYZPoint(0.,0.,0.);
00069 
00070   size_t ninputs = inputs_.size();
00071   for(size_t input = 0; input < ninputs; ++input)
00072   {
00073     edm::Handle< reco::CandidateView > h;
00074     e.getByLabel(inputs_[input], h);
00075 
00076     size_t n = h->size();
00077     for(size_t i = 0; i < n; ++i )
00078     {
00079       const reco::Candidate & c = (*h)[i];
00080 
00081       // different input collections can have different dPhi and maxZ
00082       Region r(c.momentum(), dPhi_[input], maxZ_[input]);
00083       addRegion(r);
00084     }
00085   }
00086 }
00087 
00088 
00089 void PixelUnpackingRegions::initialize(const edm::EventSetup& es)
00090 {
00091   // initialize cabling map or update it if necessary
00092   // and re-cache modules information
00093   if (watcherSiPixelFedCablingMap_.check( es ))
00094   {
00095     edm::ESTransientHandle<SiPixelFedCablingMap> cablingMap;
00096     es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
00097     cabling_.reset((SiPixelFedCabling*)cablingMap->cablingTree());
00098 
00099     edm::ESHandle<TrackerGeometry> geom;
00100     // get the TrackerGeom
00101     es.get<TrackerDigiGeometryRecord>().get( geom );
00102 
00103     phiBPIX_.clear();
00104     phiFPIXp_.clear();
00105     phiFPIXm_.clear();
00106 
00107     phiBPIX_.reserve(1024);
00108     phiFPIXp_.reserve(512);
00109     phiFPIXm_.reserve(512);
00110 
00111     std::vector<GeomDet*>::const_iterator it = geom->dets().begin();
00112     for ( ; it != geom->dets().end(); ++it)
00113     {
00114       int subdet = (*it)->geographicalId().subdetId();
00115       if (! (subdet == PixelSubdetector::PixelBarrel ||
00116              subdet == PixelSubdetector::PixelEndcap) ) continue;
00117 
00118       Module m;
00119 
00120       m.x = (*it)->position().x();
00121       m.y = (*it)->position().y();
00122       m.z = (*it)->position().z();
00123 
00124       m.phi = normalizedPhi( (*it)->position().phi() ); // ensure [-pi,+pi]
00125 
00126       m.id = (*it)->geographicalId().rawId();
00127       const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = cabling_->pathToDetUnit(m.id);
00128 
00129       m.fed = path2det[0].fed;
00130       assert(m.fed<40);
00131 
00132       if (subdet == PixelSubdetector::PixelBarrel)
00133       {
00134         phiBPIX_.push_back(m);
00135       }
00136       else if (subdet == PixelSubdetector::PixelEndcap)
00137       {
00138         if (m.z > 0.) phiFPIXp_.push_back(m);
00139         else phiFPIXm_.push_back(m);
00140       }
00141     }
00142 
00143     // pre-sort by phi
00144     std::sort(phiBPIX_.begin(),  phiBPIX_.end());
00145     std::sort(phiFPIXp_.begin(), phiFPIXp_.end());
00146     std::sort(phiFPIXm_.begin(), phiFPIXm_.end());
00147   }
00148 }
00149 
00150 
00151 void PixelUnpackingRegions::addRegion(Region &r)
00152 {
00153   ++nreg_;
00154 
00155   float phi = normalizedPhi(r.v.phi());  // ensure [-pi,+pi]
00156 
00157   Module lo(phi - r.dPhi);
00158   Module hi(phi + r.dPhi);
00159 
00160   addRegionLocal(r, phiBPIX_, lo, hi);
00161   if (r.v.eta() >  1.)
00162   {
00163     addRegionLocal(r, phiFPIXp_, lo, hi);
00164   }
00165   if (r.v.eta() < -1.)
00166   {
00167     addRegionLocal(r, phiFPIXm_, lo, hi);
00168   }
00169 }
00170 
00171 
00172 void PixelUnpackingRegions::addRegionLocal(Region &r, std::vector<Module> &container, Module lo, Module hi)
00173 {
00174   Module pi_m(-M_PI);
00175   Module pi_p( M_PI);
00176 
00177   std::vector<Module>::const_iterator a, b;
00178 
00179   if (lo.phi >= -M_PI && hi.phi <= M_PI) // interval doesn't cross the +-pi overlap
00180   {
00181     a = lower_bound(container.begin(), container.end(), lo);
00182     b = upper_bound(container.begin(), container.end(), hi);
00183     gatherFromRange(r, a, b);
00184   }
00185   else // interval is torn by the +-pi overlap
00186   {
00187     if (hi.phi >  M_PI) hi.phi -= 2.*M_PI;
00188     a = lower_bound(container.begin(), container.end(), pi_m);
00189     b = upper_bound(container.begin(), container.end(), hi);
00190     gatherFromRange(r, a, b);
00191 
00192     if (lo.phi < -M_PI) lo.phi += 2.*M_PI;
00193     a = lower_bound(container.begin(), container.end(), lo);
00194     b = upper_bound(container.begin(), container.end(), pi_p);
00195     gatherFromRange(r, a, b);
00196   }
00197 }
00198 
00199 
00200 void PixelUnpackingRegions::gatherFromRange(Region &r, std::vector<Module>::const_iterator a, std::vector<Module>::const_iterator b)
00201 {
00202   for(; a != b; ++a)
00203   {
00204     // projection in r's direction onto beam's z
00205     float zmodule = a->z - (  (a->x - beamSpot_.x())*r.cosphi + (a->y - beamSpot_.y())*r.sinphi ) * r.atantheta;
00206 
00207     // do not include modules that project too far in z
00208     if ( std::abs(zmodule) > r.maxZ ) continue;
00209 
00210     feds_.insert(a->fed);
00211     modules_.insert(a->id);
00212   }
00213 }
00214 
00215 
00216 bool PixelUnpackingRegions::mayUnpackFED(unsigned int fed_n) const
00217 {
00218   if (feds_.count(fed_n)) return true;
00219   return false;
00220 }
00221 
00222 unsigned int PixelUnpackingRegions::nBarrelFEDs() const
00223 {
00224   return std::count_if(feds_.begin(), feds_.end(), isBPIXFED );
00225 }
00226 
00227 unsigned int PixelUnpackingRegions::nForwardFEDs() const
00228 {
00229   return std::count_if(feds_.begin(), feds_.end(), isFPIXFED );
00230 }
00231 
00232 
00233 bool PixelUnpackingRegions::mayUnpackModule(unsigned int id) const
00234 {
00235   if (modules_.count(id)) return true;
00236   return false;
00237 }
00238 
00239 unsigned int PixelUnpackingRegions::nBarrelModules() const
00240 {
00241   return std::count_if(modules_.begin(), modules_.end(), isBPIXModule );
00242 }
00243 
00244 unsigned int PixelUnpackingRegions::nForwardModules() const
00245 {
00246   return std::count_if(modules_.begin(), modules_.end(), isFPIXModule );
00247 }