#include <SiStripRecHit1D.h>
Public Types | |
typedef edm::Ref < edmNew::DetSetVector < SiStripCluster > , SiStripCluster > | ClusterRef |
typedef edm::Ref < edm::LazyGetter < SiStripCluster > , SiStripCluster, edm::FindValue< SiStripCluster > > | ClusterRegionalRef |
Public Member Functions | |
virtual SiStripRecHit1D * | clone () const |
ClusterRef const & | cluster () const |
ClusterRegionalRef const & | cluster_regional () const |
virtual void | getKfComponents (KfComponentsHolder &holder) const |
bool | hasPositionAndError () const |
virtual LocalPoint | localPosition () const |
Local position. | |
virtual LocalError | localPositionError () const |
Error on the local position. | |
virtual std::vector < TrackingRecHit * > | recHits () |
Non-const access to component RecHits (if any) | |
virtual std::vector< const TrackingRecHit * > | recHits () const |
Access to component RecHits (if any) | |
void | setClusterRef (ClusterRef const &ref) |
void | setClusterRegionalRef (ClusterRegionalRef const &ref) |
void | setSigmaPitch (double sigmap) const |
virtual bool | sharesInput (const TrackingRecHit *other, SharedInputType what) const |
double | sigmaPitch () const |
SiStripRecHit1D (const LocalPoint &, const LocalError &, const DetId &, ClusterRef const &cluster) | |
SiStripRecHit1D () | |
SiStripRecHit1D (const LocalPoint &, const LocalError &, const DetId &, ClusterRegionalRef const &cluster) | |
SiStripRecHit1D (const SiStripRecHit2D *) | |
method to facilitate the convesion from 2D to 1D hits | |
~SiStripRecHit1D () | |
Private Member Functions | |
void | throwExceptionUninitialized (const char *where) const |
Private Attributes | |
ClusterRef | cluster_ |
ClusterRegionalRef | clusterRegional_ |
LocalError | err_ |
LocalPoint | pos_ |
double | sigmaPitch_ |
cache for the matcher.... |
Definition at line 13 of file SiStripRecHit1D.h.
Definition at line 21 of file SiStripRecHit1D.h.
typedef edm::Ref< edm::LazyGetter<SiStripCluster>, SiStripCluster, edm::FindValue<SiStripCluster> > SiStripRecHit1D::ClusterRegionalRef |
Definition at line 26 of file SiStripRecHit1D.h.
SiStripRecHit1D::SiStripRecHit1D | ( | ) | [inline] |
Definition at line 16 of file SiStripRecHit1D.h.
Referenced by clone().
: RecHit1D(),cluster_(),clusterRegional_(), sigmaPitch_(-1.){}
SiStripRecHit1D::~SiStripRecHit1D | ( | ) | [inline] |
Definition at line 19 of file SiStripRecHit1D.h.
{}
SiStripRecHit1D::SiStripRecHit1D | ( | const LocalPoint & | pos, |
const LocalError & | err, | ||
const DetId & | id, | ||
ClusterRef const & | cluster | ||
) |
Definition at line 5 of file SiStripRecHit1D.cc.
: RecHit1D(id),pos_(pos),err_(err), cluster_(cluster), clusterRegional_(), sigmaPitch_(-1.) {}
SiStripRecHit1D::SiStripRecHit1D | ( | const LocalPoint & | pos, |
const LocalError & | err, | ||
const DetId & | id, | ||
ClusterRegionalRef const & | cluster | ||
) |
Definition at line 15 of file SiStripRecHit1D.cc.
: RecHit1D(id),pos_(pos),err_(err), cluster_(), clusterRegional_(cluster), sigmaPitch_(-1.) {}
SiStripRecHit1D::SiStripRecHit1D | ( | const SiStripRecHit2D * | hit2D | ) |
method to facilitate the convesion from 2D to 1D hits
Definition at line 23 of file SiStripRecHit1D.cc.
References SiStripRecHit2D::cluster(), cluster_, SiStripRecHit2D::cluster_regional(), clusterRegional_, err_, edm::Ref< C, T, F >::isNonnull(), BaseSiTrackerRecHit2DLocalPos::localPositionError(), and LocalError::xx().
: RecHit1D(hit2D->geographicalId()),pos_(hit2D->localPosition()),err_(), cluster_(),clusterRegional_(),sigmaPitch_(-1) { err_ = LocalError (hit2D->localPositionError().xx(),0.,DBL_MAX); if(hit2D->cluster().isNonnull()) cluster_ = hit2D->cluster(); if(hit2D->cluster_regional().isNonnull()) clusterRegional_ = hit2D->cluster_regional(); }
virtual SiStripRecHit1D* SiStripRecHit1D::clone | ( | void | ) | const [inline, virtual] |
Implements TrackingRecHit.
Definition at line 34 of file SiStripRecHit1D.h.
References SiStripRecHit1D().
{return new SiStripRecHit1D( * this); }
ClusterRef const& SiStripRecHit1D::cluster | ( | ) | const [inline] |
Definition at line 38 of file SiStripRecHit1D.h.
References cluster_.
Referenced by fireworks::addSiStripClusters(), DeDxDiscriminatorLearner::algoAnalyze(), SiStripGainFromData::algoAnalyze(), AlignmentStats::analyze(), SiStripTrackingRecHitsValid::analyze(), TrackerHitAssociator::associateSiStripRecHit1D(), AlignmentTrackSelector::checkPrescaledHits(), reco::modules::TrackerTrackHitFilter::checkStoN(), TSiStripRecHit1D::clone(), helper::ClusterStorer::ClusterHitRecord< SiStripRecHit2D::ClusterRef >rekey(), TrackerDpgAnalysis::insertMeasurement(), AlignmentTrackSelector::isOkChargeStripHit(), TrackClusterRemover::process(), AlignmentPrescaler::produce(), TkAlCaOverlapTagger::produce(), ShallowTrackClustersProducer::produce(), ShallowGainCalibration::produce(), DeDxDiscriminatorProducer::produce(), ClusterRemovalRefSetter::reKey(), HIPAlignmentAlgorithm::run(), SiStripRecHit2D::sharesInput(), sharesInput(), and TSiStripRecHit1D::TSiStripRecHit1D().
{ return cluster_;}
ClusterRegionalRef const& SiStripRecHit1D::cluster_regional | ( | ) | const [inline] |
Definition at line 36 of file SiStripRecHit1D.h.
References clusterRegional_.
Referenced by SiStripMonitorMuonHLT::analyzeOnTrackClusters(), TrackerHitAssociator::associateSiStripRecHit1D(), TSiStripRecHit1D::clone(), SiStripRecHit2D::sharesInput(), sharesInput(), and TSiStripRecHit1D::TSiStripRecHit1D().
{ return clusterRegional_;}
void SiStripRecHit1D::getKfComponents | ( | KfComponentsHolder & | holder | ) | const [virtual] |
Reimplemented from TrackingRecHit.
Definition at line 49 of file SiStripRecHit1D.cc.
References err_, KfComponentsHolder::errors(), hasPositionAndError(), KfComponentsHolder::measuredErrors(), KfComponentsHolder::measuredParams(), KfComponentsHolder::params(), pos_, KfComponentsHolder::projection(), throwExceptionUninitialized(), KfComponentsHolder::tsosLocalErrors(), KfComponentsHolder::tsosLocalParameters(), PV3DBase< T, PVType, FrameType >::x(), and LocalError::xx().
Referenced by HelpertRecHit2DLocalPos::getKfComponents().
{ if (!hasPositionAndError()) throwExceptionUninitialized("getKfComponents"); AlgebraicVector1 & pars = holder.params<1>(); pars[0] = pos_.x(); AlgebraicSymMatrix11 & errs = holder.errors<1>(); errs(0,0) = err_.xx(); AlgebraicMatrix15 & proj = holder.projection<1>(); proj(0,3) = 1; holder.measuredParams<1>() = AlgebraicVector1( holder.tsosLocalParameters().At(3) ); holder.measuredErrors<1>() = holder.tsosLocalErrors().Sub<AlgebraicSymMatrix11>( 3, 3 ); }
bool SiStripRecHit1D::hasPositionAndError | ( | ) | const |
Definition at line 33 of file SiStripRecHit1D.cc.
References err_, pos_, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), LocalError::xy(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by getKfComponents(), localPosition(), localPositionError(), and TSiStripRecHit1D::TSiStripRecHit1D().
LocalPoint SiStripRecHit1D::localPosition | ( | ) | const [virtual] |
Local position.
Implements RecHit1D.
Definition at line 38 of file SiStripRecHit1D.cc.
References hasPositionAndError(), pos_, and throwExceptionUninitialized().
Referenced by TSiStripRecHit1D::localPosition().
{ if (!hasPositionAndError()) throwExceptionUninitialized("localPosition"); return pos_; }
LocalError SiStripRecHit1D::localPositionError | ( | ) | const [virtual] |
Error on the local position.
Implements RecHit1D.
Definition at line 43 of file SiStripRecHit1D.cc.
References err_, hasPositionAndError(), and throwExceptionUninitialized().
Referenced by TSiStripRecHit1D::localPositionError(), and TSiStripRecHit1D::parametersError().
{ if (!hasPositionAndError()) throwExceptionUninitialized("localPositionError"); return err_; }
std::vector< const TrackingRecHit * > SiStripRecHit1D::recHits | ( | ) | const [virtual] |
Access to component RecHits (if any)
Implements TrackingRecHit.
Definition at line 124 of file SiStripRecHit1D.cc.
Referenced by TSiStripRecHit1D::recHits().
{
std::vector<const TrackingRecHit*> nullvector;
return nullvector;
}
std::vector< TrackingRecHit * > SiStripRecHit1D::recHits | ( | ) | [virtual] |
Non-const access to component RecHits (if any)
Implements TrackingRecHit.
Definition at line 128 of file SiStripRecHit1D.cc.
{
std::vector<TrackingRecHit*> nullvector;
return nullvector;
}
void SiStripRecHit1D::setClusterRef | ( | ClusterRef const & | ref | ) | [inline] |
Definition at line 50 of file SiStripRecHit1D.h.
References cluster_.
Referenced by helper::ClusterStorer::ClusterHitRecord< SiStripRecHit2D::ClusterRef >rekey(), and ClusterRemovalRefSetter::reKey().
{ cluster_ = ref; }
void SiStripRecHit1D::setClusterRegionalRef | ( | ClusterRegionalRef const & | ref | ) | [inline] |
Definition at line 51 of file SiStripRecHit1D.h.
References clusterRegional_.
{ clusterRegional_ = ref; }
void SiStripRecHit1D::setSigmaPitch | ( | double | sigmap | ) | const [inline] |
bool SiStripRecHit1D::sharesInput | ( | const TrackingRecHit * | other, |
SharedInputType | what | ||
) | const [virtual] |
Returns true if the two TrackingRecHits are using the same input information (like Digis, Clusters, etc), false otherwise. The second argument specifies how much sharing is needed in order to return true: the value "all" means that all inputs of the two hits must be identical; the value "some" means that at least one of the inputs is in common.
Reimplemented from TrackingRecHit.
Definition at line 68 of file SiStripRecHit1D.cc.
References TrackingRecHit::all, cluster(), SiStripRecHit2D::cluster(), cluster_, cluster_regional(), SiStripRecHit2D::cluster_regional(), clusterRegional_, TrackingRecHit::geographicalId(), i, edm::Ref< C, T, F >::isNonnull(), TrackingRecHit::isValid(), DetId::kSubdetOffset, DetId::rawId(), and TrackingRecHit::recHits().
{ //here we exclude non si-strip subdetectors if( ((geographicalId().rawId()) >> (DetId::kSubdetOffset) ) != ( (other->geographicalId().rawId())>> (DetId::kSubdetOffset)) ) return false; //Protection against invalid hits if(! other->isValid()) return false; const std::type_info & otherType = typeid(*other); if (otherType == typeid(SiStripRecHit2D)) { const SiStripRecHit2D* otherCast = static_cast<const SiStripRecHit2D*>(other); // as 'null == null' is true, we can't just "or" the two equality tests: one of the two refs is always null! (gpetrucc) if (cluster_.isNonnull()) { return (cluster_ == otherCast->cluster()); } else { return (clusterRegional_ == otherCast->cluster_regional()); } } else if (otherType == typeid(SiStripRecHit1D)) { const SiStripRecHit1D* otherCast = static_cast<const SiStripRecHit1D*>(other); // as 'null == null' is true, we can't just "or" the two equality tests: one of the two refs is always null! (gpetrucc) if (cluster_.isNonnull()) { return (cluster_ == otherCast->cluster()); } else { return (clusterRegional_ == otherCast->cluster_regional()); } } else if (otherType == typeid(ProjectedSiStripRecHit2D)) { const SiStripRecHit2D* otherCast = & (static_cast<const ProjectedSiStripRecHit2D*>(other)->originalHit()); // as 'null == null' is true, we can't just "or" the two equality tests: one of the two refs is always null! (gpetrucc) if (cluster_.isNonnull()) { return (cluster_ == otherCast->cluster()); } else { return (clusterRegional_ == otherCast->cluster_regional()); } } else if ((otherType == typeid(SiStripMatchedRecHit2D)) && (what == all)) { return false; } else { // last resort, recur to 'recHits()', even if it returns a vector by value std::vector<const TrackingRecHit*> otherHits = other->recHits(); int ncomponents=otherHits.size(); if(ncomponents==0)return false; else if(ncomponents==1)return sharesInput(otherHits.front(),what); else if (ncomponents>1){ if(what == all )return false; else{ for(int i=0;i<ncomponents;i++){ if(sharesInput(otherHits[i],what))return true; } return false; } } return false; } }
double SiStripRecHit1D::sigmaPitch | ( | ) | const [inline] |
void SiStripRecHit1D::throwExceptionUninitialized | ( | const char * | where | ) | const [private] |
Definition at line 134 of file SiStripRecHit1D.cc.
References Exception.
Referenced by getKfComponents(), localPosition(), and localPositionError().
{ throw cms::Exception("SiStripRecHit1D") << "Trying to access " << where << " for a RecHit that was read from disk, but since CMSSW_2_1_X local positions are transient.\n" << "If you want to get coarse position/error estimation from disk, please set: ComputeCoarseLocalPositionFromDisk = True \n " << " to the TransientTrackingRecHitBuilder you are using from RecoTracker/TransientTrackingRecHit/python/TTRHBuilders_cff.py"; }
ClusterRef SiStripRecHit1D::cluster_ [private] |
Definition at line 70 of file SiStripRecHit1D.h.
Referenced by cluster(), setClusterRef(), sharesInput(), and SiStripRecHit1D().
Definition at line 73 of file SiStripRecHit1D.h.
Referenced by cluster_regional(), setClusterRegionalRef(), sharesInput(), and SiStripRecHit1D().
LocalError SiStripRecHit1D::err_ [private] |
Definition at line 67 of file SiStripRecHit1D.h.
Referenced by getKfComponents(), hasPositionAndError(), localPositionError(), and SiStripRecHit1D().
LocalPoint SiStripRecHit1D::pos_ [private] |
Definition at line 66 of file SiStripRecHit1D.h.
Referenced by getKfComponents(), hasPositionAndError(), and localPosition().
double SiStripRecHit1D::sigmaPitch_ [mutable, private] |
cache for the matcher....
Definition at line 76 of file SiStripRecHit1D.h.
Referenced by setSigmaPitch(), and sigmaPitch().