CMS 3D CMS Logo

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

PixelUnpackingRegions Class Reference

#include <PixelUnpackingRegions.h>

List of all members.

Classes

struct  Module
struct  Region

Public Member Functions

bool mayUnpackFED (unsigned int fed_n) const
 check whether a FED has to be unpacked
bool mayUnpackModule (unsigned int id) const
 check whether a module has to be unpacked
const std::set< unsigned int > * modulesToUnpack () const
 full set of module ids to unpack
unsigned int nBarrelFEDs () const
unsigned int nBarrelModules () const
unsigned int nFEDs () const
 various informational accessors:
unsigned int nForwardFEDs () const
unsigned int nForwardModules () const
unsigned int nModules () const
unsigned int nRegions () const
 PixelUnpackingRegions (const edm::ParameterSet &)
void run (const edm::Event &e, const edm::EventSetup &es)
 has to be run during each event
 ~PixelUnpackingRegions ()

Private Member Functions

void addRegion (Region &r)
void addRegionLocal (Region &r, std::vector< Module > &container, Module lo, Module hi)
void gatherFromRange (Region &r, std::vector< Module >::const_iterator, std::vector< Module >::const_iterator)
void initialize (const edm::EventSetup &es)
 run by the run method: (re)initialize the cabling data when it's necessary

Private Attributes

math::XYZPoint beamSpot_
edm::InputTag beamSpotTag_
boost::scoped_ptr
< SiPixelFedCabling
cabling_
std::vector< double > dPhi_
std::set< unsigned int > feds_
std::vector< edm::InputTaginputs_
std::vector< double > maxZ_
std::set< unsigned int > modules_
unsigned int nreg_
std::vector< ModulephiBPIX_
std::vector< ModulephiFPIXm_
std::vector< ModulephiFPIXp_
edm::ESWatcher
< SiPixelFedCablingMapRcd
watcherSiPixelFedCablingMap_

Detailed Description

Input: One or several collections of Candidate-based seeds with their objects defining the directions of unpacking regions; separate deltaPhi and maxZ tolerances could be given to each input collection. Output: FED ids and module detIds that need to be unpacked

Definition at line 27 of file PixelUnpackingRegions.h.


Constructor & Destructor Documentation

PixelUnpackingRegions::PixelUnpackingRegions ( const edm::ParameterSet conf)

Definition at line 41 of file PixelUnpackingRegions.cc.

References beamSpotTag_, dPhi_, edm::ParameterSet::getParameter(), inputs_, and maxZ_.

{
  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("Regions");
  beamSpotTag_ = regPSet.getParameter<edm::InputTag>("beamSpot");
  inputs_      = regPSet.getParameter<std::vector<edm::InputTag> >("inputs");
  dPhi_ = regPSet.getParameter<std::vector<double> >("deltaPhi");
  maxZ_ = regPSet.getParameter<std::vector<double> >("maxZ");

  if (inputs_.size() != dPhi_.size() || dPhi_.size() != maxZ_.size() )
  {
    edm::LogError("PixelUnpackingRegions")<<"Not the same size of config parameters vectors!\n"
        <<"   inputs "<<inputs_.size()<<"  deltaPhi "<<dPhi_.size() <<"  maxZ "<< maxZ_.size();
  }
}
PixelUnpackingRegions::~PixelUnpackingRegions ( ) [inline]

Definition at line 52 of file PixelUnpackingRegions.h.

{}

Member Function Documentation

void PixelUnpackingRegions::addRegion ( Region r) [private]

Definition at line 151 of file PixelUnpackingRegions.cc.

References addRegionLocal(), PixelUnpackingRegions::Region::dPhi, normalizedPhi(), nreg_, phi, phiBPIX_, phiFPIXm_, phiFPIXp_, and PixelUnpackingRegions::Region::v.

Referenced by run().

{
  ++nreg_;

  float phi = normalizedPhi(r.v.phi());  // ensure [-pi,+pi]

  Module lo(phi - r.dPhi);
  Module hi(phi + r.dPhi);

  addRegionLocal(r, phiBPIX_, lo, hi);
  if (r.v.eta() >  1.)
  {
    addRegionLocal(r, phiFPIXp_, lo, hi);
  }
  if (r.v.eta() < -1.)
  {
    addRegionLocal(r, phiFPIXm_, lo, hi);
  }
}
void PixelUnpackingRegions::addRegionLocal ( Region r,
std::vector< Module > &  container,
Module  lo,
Module  hi 
) [private]

Definition at line 172 of file PixelUnpackingRegions.cc.

References a, b, gatherFromRange(), M_PI, and PixelUnpackingRegions::Module::phi.

Referenced by addRegion().

{
  Module pi_m(-M_PI);
  Module pi_p( M_PI);

  std::vector<Module>::const_iterator a, b;

  if (lo.phi >= -M_PI && hi.phi <= M_PI) // interval doesn't cross the +-pi overlap
  {
    a = lower_bound(container.begin(), container.end(), lo);
    b = upper_bound(container.begin(), container.end(), hi);
    gatherFromRange(r, a, b);
  }
  else // interval is torn by the +-pi overlap
  {
    if (hi.phi >  M_PI) hi.phi -= 2.*M_PI;
    a = lower_bound(container.begin(), container.end(), pi_m);
    b = upper_bound(container.begin(), container.end(), hi);
    gatherFromRange(r, a, b);

    if (lo.phi < -M_PI) lo.phi += 2.*M_PI;
    a = lower_bound(container.begin(), container.end(), lo);
    b = upper_bound(container.begin(), container.end(), pi_p);
    gatherFromRange(r, a, b);
  }
}
void PixelUnpackingRegions::gatherFromRange ( Region r,
std::vector< Module >::const_iterator  a,
std::vector< Module >::const_iterator  b 
) [private]

Definition at line 200 of file PixelUnpackingRegions.cc.

References a, abs, PixelUnpackingRegions::Region::atantheta, b, beamSpot_, PixelUnpackingRegions::Region::cosphi, feds_, PixelUnpackingRegions::Region::maxZ, modules_, and PixelUnpackingRegions::Region::sinphi.

Referenced by addRegionLocal().

{
  for(; a != b; ++a)
  {
    // projection in r's direction onto beam's z
    float zmodule = a->z - (  (a->x - beamSpot_.x())*r.cosphi + (a->y - beamSpot_.y())*r.sinphi ) * r.atantheta;

    // do not include modules that project too far in z
    if ( std::abs(zmodule) > r.maxZ ) continue;

    feds_.insert(a->fed);
    modules_.insert(a->id);
  }
}
void PixelUnpackingRegions::initialize ( const edm::EventSetup es) [private]

run by the run method: (re)initialize the cabling data when it's necessary

Definition at line 89 of file PixelUnpackingRegions.cc.

References cabling_, edm::ESWatcher< T >::check(), PixelUnpackingRegions::Module::fed, relativeConstraints::geom, edm::EventSetup::get(), PixelUnpackingRegions::Module::id, m, normalizedPhi(), PixelUnpackingRegions::Module::phi, phiBPIX_, phiFPIXm_, phiFPIXp_, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, python::multivaluedict::sort(), watcherSiPixelFedCablingMap_, PixelUnpackingRegions::Module::x, PixelUnpackingRegions::Module::y, and PixelUnpackingRegions::Module::z.

Referenced by run().

{
  // initialize cabling map or update it if necessary
  // and re-cache modules information
  if (watcherSiPixelFedCablingMap_.check( es ))
  {
    edm::ESTransientHandle<SiPixelFedCablingMap> cablingMap;
    es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
    cabling_.reset((SiPixelFedCabling*)cablingMap->cablingTree());

    edm::ESHandle<TrackerGeometry> geom;
    // get the TrackerGeom
    es.get<TrackerDigiGeometryRecord>().get( geom );

    phiBPIX_.clear();
    phiFPIXp_.clear();
    phiFPIXm_.clear();

    phiBPIX_.reserve(1024);
    phiFPIXp_.reserve(512);
    phiFPIXm_.reserve(512);

    std::vector<GeomDet*>::const_iterator it = geom->dets().begin();
    for ( ; it != geom->dets().end(); ++it)
    {
      int subdet = (*it)->geographicalId().subdetId();
      if (! (subdet == PixelSubdetector::PixelBarrel ||
             subdet == PixelSubdetector::PixelEndcap) ) continue;

      Module m;

      m.x = (*it)->position().x();
      m.y = (*it)->position().y();
      m.z = (*it)->position().z();

      m.phi = normalizedPhi( (*it)->position().phi() ); // ensure [-pi,+pi]

      m.id = (*it)->geographicalId().rawId();
      const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = cabling_->pathToDetUnit(m.id);

      m.fed = path2det[0].fed;
      assert(m.fed<40);

      if (subdet == PixelSubdetector::PixelBarrel)
      {
        phiBPIX_.push_back(m);
      }
      else if (subdet == PixelSubdetector::PixelEndcap)
      {
        if (m.z > 0.) phiFPIXp_.push_back(m);
        else phiFPIXm_.push_back(m);
      }
    }

    // pre-sort by phi
    std::sort(phiBPIX_.begin(),  phiBPIX_.end());
    std::sort(phiFPIXp_.begin(), phiFPIXp_.end());
    std::sort(phiFPIXm_.begin(), phiFPIXm_.end());
  }
}
bool PixelUnpackingRegions::mayUnpackFED ( unsigned int  fed_n) const

check whether a FED has to be unpacked

Definition at line 216 of file PixelUnpackingRegions.cc.

References feds_.

Referenced by SiPixelRawToDigi::produce().

{
  if (feds_.count(fed_n)) return true;
  return false;
}
bool PixelUnpackingRegions::mayUnpackModule ( unsigned int  id) const

check whether a module has to be unpacked

Definition at line 233 of file PixelUnpackingRegions.cc.

References modules_.

{
  if (modules_.count(id)) return true;
  return false;
}
const std::set<unsigned int>* PixelUnpackingRegions::modulesToUnpack ( ) const [inline]

full set of module ids to unpack

Definition at line 64 of file PixelUnpackingRegions.h.

References modules_.

Referenced by SiPixelRawToDigi::produce().

{return &modules_;}
unsigned int PixelUnpackingRegions::nBarrelFEDs ( ) const

Definition at line 222 of file PixelUnpackingRegions.cc.

References feds_.

Referenced by SiPixelRawToDigi::produce().

{
  return std::count_if(feds_.begin(), feds_.end(), isBPIXFED );
}
unsigned int PixelUnpackingRegions::nBarrelModules ( ) const

Definition at line 239 of file PixelUnpackingRegions.cc.

References modules_.

Referenced by SiPixelRawToDigi::produce().

{
  return std::count_if(modules_.begin(), modules_.end(), isBPIXModule );
}
unsigned int PixelUnpackingRegions::nFEDs ( ) const [inline]

various informational accessors:

Definition at line 67 of file PixelUnpackingRegions.h.

References feds_.

Referenced by SiPixelRawToDigi::produce().

{ return feds_.size(); }
unsigned int PixelUnpackingRegions::nForwardFEDs ( ) const

Definition at line 227 of file PixelUnpackingRegions.cc.

References feds_.

Referenced by SiPixelRawToDigi::produce().

{
  return std::count_if(feds_.begin(), feds_.end(), isFPIXFED );
}
unsigned int PixelUnpackingRegions::nForwardModules ( ) const

Definition at line 244 of file PixelUnpackingRegions.cc.

References modules_.

Referenced by SiPixelRawToDigi::produce().

{
  return std::count_if(modules_.begin(), modules_.end(), isFPIXModule );
}
unsigned int PixelUnpackingRegions::nModules ( ) const [inline]

Definition at line 70 of file PixelUnpackingRegions.h.

References modules_.

Referenced by SiPixelRawToDigi::produce().

{ return modules_.size(); }
unsigned int PixelUnpackingRegions::nRegions ( ) const [inline]

Definition at line 73 of file PixelUnpackingRegions.h.

References nreg_.

{ return nreg_; }
void PixelUnpackingRegions::run ( const edm::Event e,
const edm::EventSetup es 
)

has to be run during each event

Definition at line 57 of file PixelUnpackingRegions.cc.

References addRegion(), SiPixelRawToDigiRegional_cfi::beamSpot, beamSpot_, beamSpotTag_, trackerHits::c, dPhi_, feds_, edm::Event::getByLabel(), h, i, initialize(), LaserDQM_cfg::input, inputs_, maxZ_, modules_, reco::Candidate::momentum(), n, nreg_, and alignCSCRings::r.

Referenced by SiPixelRawToDigi::produce().

{
  feds_.clear();
  modules_.clear();
  nreg_ = 0;

  initialize(es);

  edm::Handle<reco::BeamSpot> beamSpot;
  e.getByLabel(beamSpotTag_, beamSpot);
  beamSpot_ = beamSpot->position();
  //beamSpot_ = math::XYZPoint(0.,0.,0.);

  size_t ninputs = inputs_.size();
  for(size_t input = 0; input < ninputs; ++input)
  {
    edm::Handle< reco::CandidateView > h;
    e.getByLabel(inputs_[input], h);

    size_t n = h->size();
    for(size_t i = 0; i < n; ++i )
    {
      const reco::Candidate & c = (*h)[i];

      // different input collections can have different dPhi and maxZ
      Region r(c.momentum(), dPhi_[input], maxZ_[input]);
      addRegion(r);
    }
  }
}

Member Data Documentation

Definition at line 123 of file PixelUnpackingRegions.h.

Referenced by gatherFromRange(), and run().

Definition at line 99 of file PixelUnpackingRegions.h.

Referenced by PixelUnpackingRegions(), and run().

boost::scoped_ptr<SiPixelFedCabling> PixelUnpackingRegions::cabling_ [private]

Definition at line 122 of file PixelUnpackingRegions.h.

Referenced by initialize().

std::vector<double> PixelUnpackingRegions::dPhi_ [private]

Definition at line 97 of file PixelUnpackingRegions.h.

Referenced by PixelUnpackingRegions(), and run().

std::set<unsigned int> PixelUnpackingRegions::feds_ [private]

Definition at line 96 of file PixelUnpackingRegions.h.

Referenced by PixelUnpackingRegions(), and run().

std::vector<double> PixelUnpackingRegions::maxZ_ [private]

Definition at line 98 of file PixelUnpackingRegions.h.

Referenced by PixelUnpackingRegions(), and run().

std::set<unsigned int> PixelUnpackingRegions::modules_ [private]
unsigned int PixelUnpackingRegions::nreg_ [private]

Definition at line 103 of file PixelUnpackingRegions.h.

Referenced by addRegion(), nRegions(), and run().

std::vector<Module> PixelUnpackingRegions::phiBPIX_ [private]

Definition at line 118 of file PixelUnpackingRegions.h.

Referenced by addRegion(), and initialize().

std::vector<Module> PixelUnpackingRegions::phiFPIXm_ [private]

Definition at line 120 of file PixelUnpackingRegions.h.

Referenced by addRegion(), and initialize().

std::vector<Module> PixelUnpackingRegions::phiFPIXp_ [private]

Definition at line 119 of file PixelUnpackingRegions.h.

Referenced by addRegion(), and initialize().

Definition at line 125 of file PixelUnpackingRegions.h.

Referenced by initialize().