#include <DetIdAssociator.h>
Classes | |
struct | MapRange |
Public Types | |
typedef std::vector < GlobalPoint > ::const_iterator | const_iterator |
enum | PropagationTarget { Barrel, ForwardEndcap, BackwardEndcap } |
Public Member Functions | |
virtual void | buildMap () |
make the look-up map | |
DetIdAssociator (const int nPhi, const int nEta, const double etaBinSize) | |
double | etaBinSize () const |
look-up map bin size in eta dimension | |
virtual std::vector< DetId > | getCrossedDetIds (const std::set< DetId > &, const std::vector< GlobalPoint > &trajectory) const |
virtual std::vector< DetId > | getCrossedDetIds (const std::set< DetId > &, const std::vector< SteppingHelixStateInfo > &trajectory, const double toleranceInSigmas=-1) const |
virtual std::set< DetId > | getDetIdsCloseToAPoint (const GlobalPoint &direction, const unsigned int iNEtaPlus, const unsigned int iNEtaMinus, const unsigned int iNPhiPlus, const unsigned int iNPhiMinus) const |
virtual std::set< DetId > | getDetIdsCloseToAPoint (const GlobalPoint &direction, const MapRange &mapRange) const |
virtual std::set< DetId > | getDetIdsCloseToAPoint (const GlobalPoint &, const int iN=0) const |
virtual std::set< DetId > | getDetIdsCloseToAPoint (const GlobalPoint &point, const double d=0) const |
virtual std::set< DetId > | getDetIdsCloseToAPoint (const GlobalPoint &point, const double dThetaPlus, const double dThetaMinus, const double dPhiPlus, const double dPhiMinus) const |
virtual std::set< DetId > | getDetIdsInACone (const std::set< DetId > &, const std::vector< GlobalPoint > &trajectory, const double dR) const |
virtual const GeomDet * | getGeomDet (const DetId &) const =0 |
virtual int | iEta (const GlobalPoint &) const |
look-up map eta index | |
virtual int | iPhi (const GlobalPoint &) const |
look-up map phi index | |
virtual const char * | name () const =0 |
int | nEtaBins () const |
number of bins of the look-up map in eta dimension | |
int | nPhiBins () const |
number of bins of the look-up map in phi dimension | |
virtual void | setConditions (const DetIdAssociatorRecord &) |
virtual void | setGeometry (const DetIdAssociatorRecord &)=0 |
virtual void | setPropagator (Propagator *ptr) |
set a specific track propagator to be used | |
const FiducialVolume & | volume () const |
get active detector volume | |
virtual | ~DetIdAssociator () |
Protected Member Functions | |
virtual void | check_setup () const |
virtual bool | crossedElement (const GlobalPoint &, const GlobalPoint &, const DetId &, const double toleranceInSigmas=-1, const SteppingHelixStateInfo *=0) const |
virtual void | dumpMapContent (int, int) const |
virtual void | dumpMapContent (int, int, int, int) const |
void | fillSet (std::set< DetId > &set, unsigned int iEta, unsigned int iPhi) const |
virtual std::pair < const_iterator, const_iterator > | getDetIdPoints (const DetId &) const =0 |
virtual const unsigned int | getNumberOfSubdetectors () const |
virtual GlobalPoint | getPosition (const DetId &) const =0 |
virtual const std::vector < DetId > & | getValidDetIds (unsigned int subDetectorIndex) const =0 |
unsigned int | index (unsigned int iEta, unsigned int iPhi) const |
virtual bool | insideElement (const GlobalPoint &, const DetId &) const =0 |
virtual bool | nearElement (const GlobalPoint &point, const DetId &id, const double distance) const |
Protected Attributes | |
std::vector< DetId > | container_ |
const double | etaBinSize_ |
Propagator * | ivProp_ |
std::vector< std::pair < unsigned int, unsigned int > > | lookupMap_ |
double | maxEta_ |
double | minTheta_ |
const int | nEta_ |
const int | nPhi_ |
bool | theMapIsValid_ |
FiducialVolume | volume_ |
\
Description: Abstract base class for 3D point -> std::set<DetId>
Implementation: A look up map of active detector elements in eta-phi space is built to speed up access to the detector element geometry as well as associated hits. The map is uniformly binned in eta and phi dimensions, which can be viewed as a histogram with every bin containing DetId of elements crossed in a given eta-phi window. It is very likely that a single DetId can be found in a few bins if it's geometrical size is bigger than eta-phi bin size.
The map is implemented as a double array. The first one has fixed size and points to the range of array elements in the second one.
Definition at line 47 of file DetIdAssociator.h.
typedef std::vector<GlobalPoint>::const_iterator DetIdAssociator::const_iterator |
Definition at line 56 of file DetIdAssociator.h.
Definition at line 49 of file DetIdAssociator.h.
{ Barrel, ForwardEndcap, BackwardEndcap };
DetIdAssociator::DetIdAssociator | ( | const int | nPhi, |
const int | nEta, | ||
const double | etaBinSize | ||
) |
Definition at line 25 of file DetIdAssociator.cc.
References etaBinSize_, Exception, funct::exp(), maxEta_, minTheta_, nEta_, and nPhi_.
:nPhi_(nPhi),nEta_(nEta), lookupMap_(nPhi_*nEta_,std::pair<unsigned int, unsigned int>(0,0)), theMapIsValid_(false),etaBinSize_(etaBinSize),ivProp_(0) { if (nEta_ <= 0 || nPhi_ <= 0) throw cms::Exception("FatalError") << "incorrect look-up map size. Cannot initialize such a map."; maxEta_ = etaBinSize_*nEta_/2; minTheta_ = 2*atan(exp(-maxEta_)); }
virtual DetIdAssociator::~DetIdAssociator | ( | ) | [inline, virtual] |
Definition at line 59 of file DetIdAssociator.h.
{}
void DetIdAssociator::buildMap | ( | ) | [virtual] |
make the look-up map
Definition at line 142 of file DetIdAssociator.cc.
References abs, FiducialVolume::addActivePoint(), newFWLiteAna::bin, check_setup(), container_, FiducialVolume::determinInnerDimensions(), etaBinSize_, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, Exception, getDetIdPoints(), getNumberOfSubdetectors(), getPosition(), getValidDetIds(), iEta(), index(), info, iPhi(), edm::detail::isnan(), LogTrace, lookupMap_, FiducialVolume::maxR(), FiducialVolume::maxZ(), FiducialVolume::minR(), FiducialVolume::minZ(), name(), nEta_, nPhi_, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, launcher::step, theMapIsValid_, and volume_.
{ // the map is built in two steps to create a clean memory footprint // 0) determine how many elements each bin has // 1) fill the container map check_setup(); LogTrace("TrackAssociator")<<"building map for " << name() << "\n"; // clear the map if (nEta_ <= 0 || nPhi_ <= 0) throw cms::Exception("FatalError") << "incorrect look-up map size. Cannot build such a map."; unsigned int numberOfSubDetectors = getNumberOfSubdetectors(); unsigned totalNumberOfElementsInTheContainer(0); for ( unsigned int step = 0; step < 2; ++step){ unsigned int numberOfDetIdsOutsideEtaRange = 0; unsigned int numberOfDetIdsActive = 0; LogTrace("TrackAssociator")<< "Step: " << step; for ( unsigned int subDetectorIndex = 0; subDetectorIndex < numberOfSubDetectors; ++subDetectorIndex ){ const std::vector<DetId>& validIds = getValidDetIds(subDetectorIndex); LogTrace("TrackAssociator")<< "Number of valid DetIds for subdetector: " << subDetectorIndex << " is " << validIds.size(); for (std::vector<DetId>::const_iterator id_itr = validIds.begin(); id_itr!=validIds.end(); id_itr++) { std::pair<const_iterator,const_iterator> points = getDetIdPoints(*id_itr); LogTrace("TrackAssociatorVerbose")<< "Found " << points.second-points.first << " global points to describe geometry of DetId: " << id_itr->rawId(); int etaMax(-1); int etaMin(-1); int phiMax(-1); int phiMin(-1); // this is a bit overkill, but it should be 100% proof (when debugged :) for(std::vector<GlobalPoint>::const_iterator iter = points.first; iter != points.second; iter++) { LogTrace("TrackAssociatorVerbose")<< "\tpoint (rho,phi,z): " << iter->perp() << ", " << iter->phi() << ", " << iter->z(); // FIX ME: this should be a fatal error if(std::isnan(iter->mag())||iter->mag()>1e5) { //Detector parts cannot be 1 km away or be NaN edm::LogWarning("TrackAssociator") << "Critical error! Bad detector unit geometry:\n\tDetId:" << id_itr->rawId() << "\t mag(): " << iter->mag() << "\n" << DetIdInfo::info( *id_itr ) << "\nSkipped the element"; continue; } volume_.addActivePoint(*iter); int ieta = iEta(*iter); int iphi = iPhi(*iter); if (ieta<0 || ieta>=nEta_) { LogTrace("TrackAssociator")<<"Out of range: DetId:" << id_itr->rawId() << "\t (ieta,iphi): " << ieta << "," << iphi << "\n" << "Point: " << *iter << "\t(eta,phi): " << (*iter).eta() << "," << (*iter).phi() << "\n center: " << getPosition(*id_itr); continue; } if ( phiMin<0 ) { // first element etaMin = ieta; etaMax = ieta; phiMin = iphi; phiMax = iphi; }else{ // check for discontinuity in phi int deltaMin = abs(phiMin -iphi); int deltaMax = abs(phiMax -iphi); // assume that no single detector element has more than 3.1416 coverage in phi if ( deltaMin > nPhi_/2 && phiMin < nPhi_/2 ) phiMin+= nPhi_; if ( deltaMax > nPhi_/2 ) { if (phiMax < nPhi_/2 ) phiMax+= nPhi_; else iphi += nPhi_; } assert (iphi>=0); if ( etaMin > ieta) etaMin = ieta; if ( etaMax < ieta) etaMax = ieta; if ( phiMin > iphi) phiMin = iphi; if ( phiMax < iphi) phiMax = iphi; } } if (etaMax<0||phiMax<0||etaMin>=nEta_||phiMin>=nPhi_) { LogTrace("TrackAssociatorVerbose")<<"Out of range or no geometry: DetId:" << id_itr->rawId() << "\n\teta (min,max): " << etaMin << "," << etaMax << "\n\tphi (min,max): " << phiMin << "," << phiMax << "\nTower id: " << id_itr->rawId() << "\n"; numberOfDetIdsOutsideEtaRange++; continue; } numberOfDetIdsActive++; LogTrace("TrackAssociatorVerbose") << "DetId (ieta_min,ieta_max,iphi_min,iphi_max): " << id_itr->rawId() << ", " << etaMin << ", " << etaMax << ", " << phiMin << ", " << phiMax; for(int ieta = etaMin; ieta <= etaMax; ieta++) for(int iphi = phiMin; iphi <= phiMax; iphi++) if ( step == 0 ){ lookupMap_.at(index(ieta,iphi%nPhi_)).second++; totalNumberOfElementsInTheContainer++; } else { container_.at(lookupMap_.at(index(ieta,iphi%nPhi_)).first) = *id_itr; lookupMap_.at(index(ieta,iphi%nPhi_)).first--; totalNumberOfElementsInTheContainer--; } } } LogTrace("TrackAssociator") << "Number of elements outside the allowed range ( |eta|>"<< nEta_/2*etaBinSize_ << "): " << numberOfDetIdsOutsideEtaRange << "\n"; LogTrace("TrackAssociator") << "Number of active DetId's mapped: " << numberOfDetIdsActive << "\n"; if ( step == 0 ){ // allocate container_.resize(totalNumberOfElementsInTheContainer); // fill the range index in the lookup map to point to the last element in the range unsigned int index(0); for ( std::vector<std::pair<unsigned int,unsigned int> >::iterator bin = lookupMap_.begin(); bin != lookupMap_.end(); ++bin ) { if (bin->second==0) continue; index += bin->second; bin->first = index-1; } } } if ( totalNumberOfElementsInTheContainer != 0 ) throw cms::Exception("FatalError") << "Look-up map filled incorrectly. Structural problem. Get in touch with the developer."; volume_.determinInnerDimensions(); edm::LogVerbatim("TrackAssociator") << "Fiducial volume for " << name() << " (minR, maxR, minZ, maxZ): " << volume_.minR() << ", " << volume_.maxR() << ", " << volume_.minZ() << ", " << volume_.maxZ(); theMapIsValid_ = true; }
void DetIdAssociator::check_setup | ( | ) | const [protected, virtual] |
Reimplemented in CaloDetIdAssociator, and MuonDetIdAssociator.
Definition at line 372 of file DetIdAssociator.cc.
References etaBinSize_, Exception, nEta_, and nPhi_.
Referenced by buildMap(), getCrossedDetIds(), getDetIdsCloseToAPoint(), and getDetIdsInACone().
{ if (nEta_==0) throw cms::Exception("FatalError") << "Number of eta bins is not set.\n"; if (nPhi_==0) throw cms::Exception("FatalError") << "Number of phi bins is not set.\n"; // if (ivProp_==0) throw cms::Exception("FatalError") << "Track propagator is not defined\n"; if (etaBinSize_==0) throw cms::Exception("FatalError") << "Eta bin size is not set.\n"; }
virtual bool DetIdAssociator::crossedElement | ( | const GlobalPoint & | , |
const GlobalPoint & | , | ||
const DetId & | , | ||
const double | toleranceInSigmas = -1 , |
||
const SteppingHelixStateInfo * | = 0 |
||
) | const [inline, protected, virtual] |
Reimplemented in CaloDetIdAssociator.
Definition at line 137 of file DetIdAssociator.h.
Referenced by getCrossedDetIds().
{ return false; }
void DetIdAssociator::dumpMapContent | ( | int | ieta, |
int | iphi | ||
) | const [protected, virtual] |
Definition at line 320 of file DetIdAssociator.cc.
References begin, fillSet(), getDetIdPoints(), LogTrace, nPhi_, and point.
Referenced by dumpMapContent().
{ if (! (ieta>=0 && ieta<nEta_ && iphi>=0) ) { edm::LogWarning("TrackAssociator") << "ieta or iphi is out of range. Skipped."; return; } std::set<DetId> set; fillSet(set,ieta,iphi%nPhi_); LogTrace("TrackAssociator") << "Map content for cell (ieta,iphi): " << ieta << ", " << iphi%nPhi_; for(std::set<DetId>::const_iterator itr = set.begin(); itr!=set.end(); itr++) { LogTrace("TrackAssociator") << "\tDetId " << itr->rawId() << ", geometry (x,y,z,rho,eta,phi):"; std::pair<const_iterator,const_iterator> points = getDetIdPoints(*itr); for(std::vector<GlobalPoint>::const_iterator point = points.first; point != points.second; point++) LogTrace("TrackAssociator") << "\t\t" << point->x() << ", " << point->y() << ", " << point->z() << ", " << point->perp() << ", " << point->eta() << ", " << point->phi(); } }
void DetIdAssociator::dumpMapContent | ( | int | ieta_min, |
int | ieta_max, | ||
int | iphi_min, | ||
int | iphi_max | ||
) | const [protected, virtual] |
Definition at line 341 of file DetIdAssociator.cc.
References dumpMapContent(), i, and j.
double DetIdAssociator::etaBinSize | ( | ) | const [inline] |
look-up map bin size in eta dimension
Definition at line 112 of file DetIdAssociator.h.
References etaBinSize_.
{ return etaBinSize_;};
void DetIdAssociator::fillSet | ( | std::set< DetId > & | set, |
unsigned int | iEta, | ||
unsigned int | iPhi | ||
) | const [protected] |
Definition at line 380 of file DetIdAssociator.cc.
References container_, i, index(), edm::eventsetup::heterocontainer::insert(), lookupMap_, and findQualityFiles::size.
Referenced by dumpMapContent(), and getDetIdsCloseToAPoint().
{ unsigned int i = index(iEta,iPhi); unsigned int i0 = lookupMap_.at(i).first; unsigned int size = lookupMap_.at(i).second; for ( i = i0; i < i0+size; ++i ) set.insert(container_.at(i)); }
std::vector< DetId > DetIdAssociator::getCrossedDetIds | ( | const std::set< DetId > & | inset, |
const std::vector< SteppingHelixStateInfo > & | trajectory, | ||
const double | toleranceInSigmas = -1 |
||
) | const [virtual] |
Definition at line 299 of file DetIdAssociator.cc.
References check_setup(), crossedElement(), i, convertSQLitetoXML_cfg::output, and position.
{ check_setup(); std::vector<DetId> output; std::set<DetId> ids(inset); for ( unsigned int i=0; i+1 < trajectory.size(); ++i ) { std::set<DetId>::const_iterator id_iter = ids.begin(); while ( id_iter != ids.end() ) { if ( crossedElement(trajectory[i].position(),trajectory[i+1].position(),*id_iter,tolerance,&trajectory[i]) ){ output.push_back(*id_iter); ids.erase(id_iter++); }else id_iter++; } } return output; }
std::vector< DetId > DetIdAssociator::getCrossedDetIds | ( | const std::set< DetId > & | inset, |
const std::vector< GlobalPoint > & | trajectory | ||
) | const [virtual] |
Definition at line 280 of file DetIdAssociator.cc.
References check_setup(), crossedElement(), i, and convertSQLitetoXML_cfg::output.
{ check_setup(); std::vector<DetId> output; std::set<DetId> ids(inset); for ( unsigned int i=0; i+1 < trajectory.size(); ++i ) { std::set<DetId>::const_iterator id_iter = ids.begin(); while ( id_iter != ids.end() ) { if ( crossedElement(trajectory[i],trajectory[i+1],*id_iter) ){ output.push_back(*id_iter); ids.erase(id_iter++); }else id_iter++; } } return output; }
virtual std::pair<const_iterator, const_iterator> DetIdAssociator::getDetIdPoints | ( | const DetId & | ) | const [protected, pure virtual] |
Implemented in CaloDetIdAssociator, and MuonDetIdAssociator.
Referenced by buildMap(), and dumpMapContent().
std::set< DetId > DetIdAssociator::getDetIdsCloseToAPoint | ( | const GlobalPoint & | direction, |
const unsigned int | iNEtaPlus, | ||
const unsigned int | iNEtaMinus, | ||
const unsigned int | iNPhiPlus, | ||
const unsigned int | iNPhiMinus | ||
) | const [virtual] |
Definition at line 43 of file DetIdAssociator.cc.
References check_setup(), PV3DBase< T, PVType, FrameType >::eta(), Exception, fillSet(), i, iEta(), iPhi(), j, LogTrace, nEta_, nPhi_, PV3DBase< T, PVType, FrameType >::phi(), and theMapIsValid_.
{ std::set<DetId> set; check_setup(); if (! theMapIsValid_ ) throw cms::Exception("FatalError") << "map is not valid."; LogTrace("TrackAssociator") << "(iNEtaPlus, iNEtaMinus, iNPhiPlus, iNPhiMinus): " << iNEtaPlus << ", " << iNEtaMinus << ", " << iNPhiPlus << ", " << iNPhiMinus; LogTrace("TrackAssociator") << "point (eta,phi): " << direction.eta() << "," << direction.phi(); int ieta = iEta(direction); int iphi = iPhi(direction); LogTrace("TrackAssociator") << "(ieta,iphi): " << ieta << "," << iphi << "\n"; if (ieta>=0 && ieta<nEta_ && iphi>=0 && iphi<nPhi_){ fillSet(set,ieta,iphi); // dumpMapContent(ieta,iphi); // check if any neighbor bin is requested if (iNEtaPlus + iNEtaMinus + iNPhiPlus + iNPhiMinus >0 ){ LogTrace("TrackAssociator") << "Add neighbors (ieta,iphi): " << ieta << "," << iphi; // eta int maxIEta = ieta+iNEtaPlus; int minIEta = ieta-iNEtaMinus; if (maxIEta>=nEta_) maxIEta = nEta_-1; if (minIEta<0) minIEta = 0; // phi int maxIPhi = iphi+iNPhiPlus; int minIPhi = iphi-iNPhiMinus; if (maxIPhi-minIPhi>=nPhi_){ // all elements in phi minIPhi = 0; maxIPhi = nPhi_-1; } if(minIPhi<0) { minIPhi+=nPhi_; maxIPhi+=nPhi_; } LogTrace("TrackAssociator") << "\tieta (min,max): " << minIEta << "," << maxIEta; LogTrace("TrackAssociator") << "\tiphi (min,max): " << minIPhi << "," << maxIPhi<< "\n"; // dumpMapContent(minIEta,maxIEta,minIPhi,maxIPhi); for (int i=minIEta;i<=maxIEta;i++) for (int j=minIPhi;j<=maxIPhi;j++) { if( i==ieta && j==iphi) continue; // already in the set fillSet(set,i,j%nPhi_); } } } return set; }
std::set< DetId > DetIdAssociator::getDetIdsCloseToAPoint | ( | const GlobalPoint & | point, |
const double | d = 0 |
||
) | const [virtual] |
Preselect DetIds close to a point on the inner surface of the detector. "d" defines the allowed range in theta-phi space:
Definition at line 94 of file DetIdAssociator.cc.
References getDetIdsCloseToAPoint().
{ return getDetIdsCloseToAPoint(point,d,d,d,d); }
std::set< DetId > DetIdAssociator::getDetIdsCloseToAPoint | ( | const GlobalPoint & | point, |
const double | dThetaPlus, | ||
const double | dThetaMinus, | ||
const double | dPhiPlus, | ||
const double | dPhiMinus | ||
) | const [virtual] |
Definition at line 100 of file DetIdAssociator.cc.
References abs, PV3DBase< T, PVType, FrameType >::eta(), etaBinSize_, getDetIdsCloseToAPoint(), funct::log(), LogTrace, M_PI, maxEta, benchmark_cfg::minEta, minTheta_, n, nPhi_, funct::tan(), and PV3DBase< T, PVType, FrameType >::theta().
{ LogTrace("TrackAssociator") << "(dThetaPlus,dThetaMinus,dPhiPlus,dPhiMinus): " << dThetaPlus << ", " << dThetaMinus << ", " << dPhiPlus << ", " << dPhiMinus; unsigned int n = 0; if ( dThetaPlus<0 || dThetaMinus<0 || dPhiPlus<0 || dPhiMinus<0) return getDetIdsCloseToAPoint(point,n,n,n,n); // check that region of interest overlaps with the look-up map double maxTheta = point.theta()+dThetaPlus; if (maxTheta > M_PI-minTheta_) maxTheta = M_PI-minTheta_; double minTheta = point.theta()-dThetaMinus; if (minTheta < minTheta_) minTheta = minTheta_; if ( maxTheta < minTheta_ || minTheta > M_PI-minTheta_) return std::set<DetId>(); // take into account non-linear dependence of eta from // theta in regions with large |eta| double minEta = -log(tan(maxTheta/2)); double maxEta = -log(tan(minTheta/2)); unsigned int iNEtaPlus = abs(int( ( maxEta-point.eta() )/etaBinSize_)); unsigned int iNEtaMinus = abs(int( ( point.eta() - minEta )/etaBinSize_)); unsigned int iNPhiPlus = abs(int( dPhiPlus/(2*M_PI)*nPhi_ )); unsigned int iNPhiMinus = abs(int( dPhiMinus/(2*M_PI)*nPhi_ )); // add one more bin in each direction to guaranty that we don't miss anything return getDetIdsCloseToAPoint(point, iNEtaPlus+1, iNEtaMinus+1, iNPhiPlus+1, iNPhiMinus+1); }
std::set< DetId > DetIdAssociator::getDetIdsCloseToAPoint | ( | const GlobalPoint & | direction, |
const MapRange & | mapRange | ||
) | const [virtual] |
Definition at line 354 of file DetIdAssociator.cc.
References DetIdAssociator::MapRange::dPhiMinus, DetIdAssociator::MapRange::dPhiPlus, DetIdAssociator::MapRange::dThetaMinus, DetIdAssociator::MapRange::dThetaPlus, and getDetIdsCloseToAPoint().
{ return getDetIdsCloseToAPoint(direction, mapRange.dThetaPlus, mapRange.dThetaMinus, mapRange.dPhiPlus, mapRange.dPhiMinus); }
std::set< DetId > DetIdAssociator::getDetIdsCloseToAPoint | ( | const GlobalPoint & | direction, |
const int | iN = 0 |
||
) | const [virtual] |
Preselect DetIds close to a point on the inner surface of the detector. "iN" is a number of the adjacent bins of the map to retrieve
Definition at line 35 of file DetIdAssociator.cc.
References n.
Referenced by getDetIdsCloseToAPoint().
{ unsigned int n = 0; if (iN>0) n = iN; return getDetIdsCloseToAPoint(direction,n,n,n,n); }
std::set< DetId > DetIdAssociator::getDetIdsInACone | ( | const std::set< DetId > & | inset, |
const std::vector< GlobalPoint > & | trajectory, | ||
const double | dR | ||
) | const [virtual] |
Find DetIds that satisfy given requirements
Definition at line 264 of file DetIdAssociator.cc.
References check_setup(), M_PI, maxEta_, and nearElement().
{ if ( dR > 2*M_PI && dR > maxEta_ ) return inset; check_setup(); std::set<DetId> outset; for(std::set<DetId>::const_iterator id_iter = inset.begin(); id_iter != inset.end(); id_iter++) for(std::vector<GlobalPoint>::const_iterator point_iter = trajectory.begin(); point_iter != trajectory.end(); point_iter++) if (nearElement(*point_iter,*id_iter,dR)) { outset.insert(*id_iter); break; } return outset; }
Implemented in CaloDetIdAssociator, and MuonDetIdAssociator.
virtual const unsigned int DetIdAssociator::getNumberOfSubdetectors | ( | ) | const [inline, protected, virtual] |
Reimplemented in EcalDetIdAssociator, and HcalDetIdAssociator.
Definition at line 132 of file DetIdAssociator.h.
Referenced by buildMap().
{ return 1;}
virtual GlobalPoint DetIdAssociator::getPosition | ( | const DetId & | ) | const [protected, pure virtual] |
Implemented in CaloDetIdAssociator, and MuonDetIdAssociator.
Referenced by buildMap(), and nearElement().
virtual const std::vector<DetId>& DetIdAssociator::getValidDetIds | ( | unsigned int | subDetectorIndex | ) | const [protected, pure virtual] |
Implemented in CaloDetIdAssociator, EcalDetIdAssociator, HcalDetIdAssociator, HODetIdAssociator, MuonDetIdAssociator, and PreshowerDetIdAssociator.
Referenced by buildMap().
int DetIdAssociator::iEta | ( | const GlobalPoint & | point | ) | const [virtual] |
look-up map eta index
Definition at line 131 of file DetIdAssociator.cc.
References PV3DBase< T, PVType, FrameType >::eta(), etaBinSize_, and nEta_.
Referenced by buildMap(), and getDetIdsCloseToAPoint().
{ return int(point.eta()/etaBinSize_ + nEta_/2); }
unsigned int DetIdAssociator::index | ( | unsigned int | iEta, |
unsigned int | iPhi | ||
) | const [inline, protected] |
Definition at line 146 of file DetIdAssociator.h.
Referenced by buildMap(), and fillSet().
virtual bool DetIdAssociator::insideElement | ( | const GlobalPoint & | , |
const DetId & | |||
) | const [protected, pure virtual] |
Implemented in CaloDetIdAssociator, and MuonDetIdAssociator.
int DetIdAssociator::iPhi | ( | const GlobalPoint & | point | ) | const [virtual] |
look-up map phi index
Definition at line 136 of file DetIdAssociator.cc.
References M_PI, nPhi_, and PV3DBase< T, PVType, FrameType >::phi().
Referenced by buildMap(), getDetIdsCloseToAPoint(), and index().
virtual const char* DetIdAssociator::name | ( | ) | const [pure virtual] |
Implemented in CaloDetIdAssociator, EcalDetIdAssociator, HcalDetIdAssociator, HODetIdAssociator, MuonDetIdAssociator, and PreshowerDetIdAssociator.
Referenced by buildMap().
bool DetIdAssociator::nearElement | ( | const GlobalPoint & | point, |
const DetId & | id, | ||
const double | distance | ||
) | const [protected, virtual] |
Definition at line 362 of file DetIdAssociator.cc.
References SiPixelRawToDigiRegional_cfi::deltaPhi, PV3DBase< T, PVType, FrameType >::eta(), getPosition(), M_PI, and PV3DBase< T, PVType, FrameType >::phi().
Referenced by getDetIdsInACone().
int DetIdAssociator::nEtaBins | ( | ) | const [inline] |
number of bins of the look-up map in eta dimension
Definition at line 110 of file DetIdAssociator.h.
References nEta_.
{ return nEta_;}
int DetIdAssociator::nPhiBins | ( | ) | const [inline] |
number of bins of the look-up map in phi dimension
Definition at line 108 of file DetIdAssociator.h.
References nPhi_.
{ return nPhi_;}
virtual void DetIdAssociator::setConditions | ( | const DetIdAssociatorRecord & | ) | [inline, virtual] |
virtual void DetIdAssociator::setGeometry | ( | const DetIdAssociatorRecord & | ) | [pure virtual] |
Implemented in CaloDetIdAssociator, and MuonDetIdAssociator.
virtual void DetIdAssociator::setPropagator | ( | Propagator * | ptr | ) | [inline, virtual] |
set a specific track propagator to be used
Definition at line 106 of file DetIdAssociator.h.
References ivProp_.
{ ivProp_ = ptr; };
const FiducialVolume & DetIdAssociator::volume | ( | void | ) | const |
get active detector volume
Definition at line 348 of file DetIdAssociator.cc.
References Exception, theMapIsValid_, and volume_.
Referenced by BetaCalculatorECAL::calcEcalDeposit().
{ if (! theMapIsValid_ ) throw cms::Exception("FatalError") << "map is not valid."; return volume_; }
std::vector<DetId> DetIdAssociator::container_ [protected] |
Definition at line 158 of file DetIdAssociator.h.
Referenced by buildMap(), and fillSet().
const double DetIdAssociator::etaBinSize_ [protected] |
Definition at line 160 of file DetIdAssociator.h.
Referenced by buildMap(), check_setup(), DetIdAssociator(), etaBinSize(), getDetIdsCloseToAPoint(), and iEta().
Propagator* DetIdAssociator::ivProp_ [protected] |
Definition at line 164 of file DetIdAssociator.h.
Referenced by setPropagator().
std::vector<std::pair<unsigned int,unsigned int> > DetIdAssociator::lookupMap_ [protected] |
Definition at line 157 of file DetIdAssociator.h.
Referenced by buildMap(), and fillSet().
double DetIdAssociator::maxEta_ [protected] |
Definition at line 161 of file DetIdAssociator.h.
Referenced by DetIdAssociator(), and getDetIdsInACone().
double DetIdAssociator::minTheta_ [protected] |
Definition at line 162 of file DetIdAssociator.h.
Referenced by DetIdAssociator(), and getDetIdsCloseToAPoint().
const int DetIdAssociator::nEta_ [protected] |
Definition at line 153 of file DetIdAssociator.h.
Referenced by buildMap(), check_setup(), DetIdAssociator(), getDetIdsCloseToAPoint(), iEta(), and nEtaBins().
const int DetIdAssociator::nPhi_ [protected] |
Definition at line 152 of file DetIdAssociator.h.
Referenced by buildMap(), check_setup(), DetIdAssociator(), dumpMapContent(), getDetIdsCloseToAPoint(), index(), iPhi(), and nPhiBins().
bool DetIdAssociator::theMapIsValid_ [protected] |
Definition at line 159 of file DetIdAssociator.h.
Referenced by buildMap(), getDetIdsCloseToAPoint(), and volume().
FiducialVolume DetIdAssociator::volume_ [protected] |
Definition at line 169 of file DetIdAssociator.h.
Referenced by buildMap(), and volume().