CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
reco::modules::TrackerTrackHitFilter Class Reference
Inheritance diagram for reco::modules::TrackerTrackHitFilter:
edm::stream::EDProducer<> edm::stream::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

class  Rule
 

Public Member Functions

int checkHit (const edm::EventSetup &iSetup, const DetId &detid, const TrackingRecHit *hit)
 
virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
void produceFromTrack (const edm::EventSetup &iSetup, const Track *itt, std::vector< TrackingRecHit * > &hits)
 
void produceFromTrajectory (const edm::EventSetup &iSetup, const Trajectory *itt, std::vector< TrackingRecHit * > &hits)
 
 TrackerTrackHitFilter (const edm::ParameterSet &iConfig)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::stream::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

bool checkHitAngle (const TrajectoryMeasurement &meas)
 
bool checkPXLCorrClustCharge (const TrajectoryMeasurement &meas)
 
bool checkStoN (const edm::EventSetup &iSetup, const DetId &id, const TrackingRecHit *therechit)
 
int layerFromId (const DetId &id, const TrackerTopology *tTopo) const
 
TrackCandidate makeCandidate (const reco::Track &tk, std::vector< TrackingRecHit * >::iterator hitsBegin, std::vector< TrackingRecHit * >::iterator hitsEnd)
 
void parseStoN (const std::string &str)
 
int sideFromId (const DetId &id, const TrackerTopology *tTopo) const
 

Private Attributes

bool checkPXLQuality_
 
std::string CMNSubtractionMode_
 
std::vector< uint32_t > detsToIgnore_
 
edm::EventNumber_t iEvt
 
edm::RunNumber_t iRun
 
size_t minimumHits_
 
int nOverlaps
 
double PXLcorrClusChargeCut_
 
double pxlTPLProbXY_
 
double pxlTPLProbXYQ_
 
std::vector< int32_t > pxlTPLqBin_
 
bool rejectBadStoNHits_
 
bool rejectLowAngleHits_
 
bool replaceWithInactiveHits_
 
std::vector< Rulerules_
 
edm::InputTag src_
 
bool stripAllInvalidHits_
 
bool stripBackInvalidHits_
 
bool stripFrontInvalidHits_
 
std::vector< bool > subdetStoN_
 
std::vector< double > subdetStoNhighcut_
 
std::vector< double > subdetStoNlowcut_
 
bool tagOverlaps_
 
edm::ESHandle
< TransientTrackingRecHitBuilder
theBuilder
 
edm::ESHandle< TrackerGeometrytheGeometry
 
edm::ESHandle< MagneticFieldtheMagField
 
edm::EDGetTokenT
< reco::TrackCollection
tokenTracks
 
edm::EDGetTokenT
< TrajTrackAssociationCollection
tokenTrajTrack
 
double TrackAngleCut_
 
bool useTrajectories_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T...> CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T...> HasAbility
 
typedef
CacheTypes::LuminosityBlockCache 
LuminosityBlockCache
 
typedef
LuminosityBlockContextT
< LuminosityBlockCache,
RunCache, GlobalCache
LuminosityBlockContext
 
typedef
CacheTypes::LuminosityBlockSummaryCache 
LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache,
GlobalCache
RunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDProducerBase
typedef EDProducerAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::stream::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 75 of file TrackerTrackHitFilter.cc.

Constructor & Destructor Documentation

reco::modules::TrackerTrackHitFilter::TrackerTrackHitFilter ( const edm::ParameterSet iConfig)

Definition at line 275 of file TrackerTrackHitFilter.cc.

References detsToIgnore_, Exception, edm::ParameterSet::getParameter(), parseStoN(), PXLcorrClusChargeCut_, pxlTPLqBin_, rejectBadStoNHits_, rejectLowAngleHits_, replaceWithInactiveHits_, rules_, python.multivaluedict::sort(), src_, stripAllInvalidHits_, subdetStoN_, subdetStoNhighcut_, subdetStoNlowcut_, tokenTracks, tokenTrajTrack, TrackAngleCut_, and useTrajectories_.

275  :
276  src_(iConfig.getParameter<edm::InputTag>("src")),
277  minimumHits_(iConfig.getParameter<uint32_t>("minimumHits")),
278  replaceWithInactiveHits_(iConfig.getParameter<bool>("replaceWithInactiveHits")),
279  stripFrontInvalidHits_(iConfig.getParameter<bool>("stripFrontInvalidHits")),
280  stripBackInvalidHits_( iConfig.getParameter<bool>("stripBackInvalidHits") ),
281  stripAllInvalidHits_( iConfig.getParameter<bool>("stripAllInvalidHits") ),
282  rejectBadStoNHits_( iConfig.getParameter<bool>("rejectBadStoNHits") ),
283  CMNSubtractionMode_( iConfig.getParameter<std::string>("CMNSubtractionMode") ),
284  detsToIgnore_( iConfig.getParameter<std::vector<uint32_t> >("detsToIgnore") ),
285  useTrajectories_( iConfig.getParameter<bool>("useTrajectories") ),
286  rejectLowAngleHits_( iConfig.getParameter<bool>("rejectLowAngleHits") ),
287  TrackAngleCut_( iConfig.getParameter<double>("TrackAngleCut") ),
288  checkPXLQuality_(iConfig.getParameter<bool>("usePixelQualityFlag") ),
289  pxlTPLProbXY_(iConfig.getParameter<double>("PxlTemplateProbXYCut")),
290  pxlTPLProbXYQ_(iConfig.getParameter<double>("PxlTemplateProbXYChargeCut")),
291  pxlTPLqBin_(iConfig.getParameter<std::vector<int32_t> >("PxlTemplateqBinCut")),
292  PXLcorrClusChargeCut_(iConfig.getParameter<double>("PxlCorrClusterChargeCut")),
293  tagOverlaps_( iConfig.getParameter<bool>("tagOverlaps") )
294 {
295 
296  if(useTrajectories_) tokenTrajTrack = consumes<TrajTrackAssociationCollection>(src_);
297  else tokenTracks = consumes<reco::TrackCollection>(src_);
298 
299  // sanity check
301  throw cms::Exception("Configuration") << "Inconsistent Configuration: you can't set both 'stripAllInvalidHits' and 'replaceWithInactiveHits' to true\n";
302  }
304  throw cms::Exception("Configuration") << "Wrong configuration of TrackerTrackHitFilter. You cannot apply the cut on the track angle without using Trajectories!\n";
305  }
307  throw cms::Exception("Configuration") << "Wrong configuration of TrackerTrackHitFilter. You cannot apply the cut on the corrected pixel cluster charge without using Trajectories!\n";
308  }
309 
310  if(pxlTPLqBin_.size()>2){
311  edm::LogInfo("TrackerTrackHitFIlter")<<"Warning from TrackerTrackHitFilter: vector with qBin cuts has size > 2. Additional items will be ignored.";
312  }
313 
314 
315 
316  // read and parse commands
317  std::vector<std::string> str_rules = iConfig.getParameter<std::vector<std::string> >("commands");
318  rules_.reserve(str_rules.size());
319  for (std::vector<std::string>::const_iterator it = str_rules.begin(), ed = str_rules.end(); it != ed; ++it) {
320  rules_.push_back(Rule(*it));
321  }
322 
323  if(rejectBadStoNHits_){//commands for S/N cut
324 
325  subdetStoN_.reserve(6);
326  subdetStoNlowcut_.reserve(6);
327  subdetStoNhighcut_.reserve(6);
328  int cnt=0;
329  for(cnt=0;cnt<6;cnt++ ){
330  subdetStoN_[cnt]=false;
331  subdetStoNlowcut_[cnt]=-1.0;
332  subdetStoNhighcut_[cnt]=-1.0;
333  }
334 
335  std::vector<std::string> str_StoNrules = iConfig.getParameter<std::vector<std::string> >("StoNcommands");
336  for (std::vector<std::string>::const_iterator str_StoN = str_StoNrules.begin(); str_StoN != str_StoNrules.end(); ++str_StoN) {
337  parseStoN(*str_StoN);
338  }
340  edm::LogInfo("TrackerTrackHitFilter")<<"Finished parsing S/N. Applying following cuts to subdets:";
341  for(cnt=0;cnt<6;cnt++ ){
343  edm::LogVerbatim("TrackerTrackHitFilter")<<"Subdet #"<<cnt+1<<" -> "<<subdetStoNlowcut_[cnt]<<" , "<<subdetStoNhighcut_[cnt];
344  }
345  }//end if rejectBadStoNHits_
346 
347 
348  if(rejectLowAngleHits_ ) edm::LogInfo("TrackerTrackHitFilter")<<"\nApplying cut on angle track = "<<TrackAngleCut_;
349 
350 
351  // sort detids to ignore
352  std::sort(detsToIgnore_.begin(), detsToIgnore_.end());
353 
354  // issue the produce<>
355  produces<TrackCandidateCollection>();
356 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > tokenTracks
void parseStoN(const std::string &str)
edm::EDGetTokenT< TrajTrackAssociationCollection > tokenTrajTrack

Member Function Documentation

int reco::modules::TrackerTrackHitFilter::checkHit ( const edm::EventSetup iSetup,
const DetId detid,
const TrackingRecHit hit 
)

Definition at line 702 of file TrackerTrackHitFilter.cc.

References checkStoN(), DetId::det(), detsToIgnore_, edm::EventSetup::get(), TrackingRecHit::isValid(), edm::ESHandle< class >::product(), DetId::rawId(), rejectBadStoNHits_, rules_, and DetId::Tracker.

Referenced by produceFromTrack(), and produceFromTrajectory().

702  {
703 
704  //Retrieve tracker topology from geometry
706  iSetup.get<TrackerTopologyRcd>().get(tTopoHand);
707  const TrackerTopology *tTopo=tTopoHand.product();
708 
709  int hitresult=0;
710  if (hit->isValid()) {
711 
712  if (detid.det() == DetId::Tracker) { // check for tracker hits
713  bool verdict = true;
714  // first check at structure level
715  for (std::vector<Rule>::const_iterator itr = rules_.begin(), edr = rules_.end(); itr != edr; ++itr) {
716  itr->apply(detid, tTopo, verdict);
717  }
718 
719  // if the hit is good, check again at module level
720  if ( verdict){
721  if(std::binary_search(detsToIgnore_.begin(), detsToIgnore_.end(), detid.rawId())) {
722  hitresult=-4;
723  }
724  }
725  else hitresult=-3;
726  //if the hit is in the desired part of the det, check other things
727  if( hitresult==0 && rejectBadStoNHits_){
728  if(!checkStoN(iSetup, detid, hit))hitresult=-5;
729  }//end if S/N is ok
730  }//end hit in tracker
731  else hitresult =-2;
732  }//end hit is valid
733  else hitresult = -1; //invalid hit
734  return hitresult;
735 }//end TrackerTrackHitFilter::checkHit()
bool checkStoN(const edm::EventSetup &iSetup, const DetId &id, const TrackingRecHit *therechit)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool isValid() const
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
bool reco::modules::TrackerTrackHitFilter::checkHitAngle ( const TrajectoryMeasurement meas)
private

Definition at line 834 of file TrackerTrackHitFilter.cc.

References Abs(), angle(), SiPixelRecHit::cluster(), TrackingRecHit::geographicalId(), GeomDetEnumerators::isTrackerPixel(), TrajectoryStateOnSurface::isValid(), TrajectoryStateOnSurface::localDirection(), funct::pow(), PXLcorrClusChargeCut_, TrajectoryMeasurement::recHit(), rejectLowAngleHits_, mathSSE::sqrt(), DetId::subdetId(), funct::tan(), theGeometry, TrackAngleCut_, TrajectoryMeasurement::updatedState(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by produceFromTrajectory().

834  {
835 
836  bool angle_ok=false;
837  bool corrcharge_ok=true;
839  /*
840  edm::LogDebug("TrackerTrackHitFilter")<<"TSOS parameters: ";
841  edm::LogDebug("TrackerTrackHitFilter") <<"Global momentum: "<<tsos.globalMomentum().x()<<" "<<tsos.globalMomentum().y()<<" "<<tsos.globalMomentum().z();
842  edm::LogDebug("TrackerTrackHitFilter") <<"Local momentum: "<<tsos.localMomentum().x()<<" "<<tsos.localMomentum().y()<<" "<<tsos.localMomentum().z();
843  edm::LogDebug("TrackerTrackHitFilter") <<"Track charge: " <<tsos.charge();
844  edm::LogDebug("TrackerTrackHitFilter")<<"Local position: " <<tsos.localPosition().x()<<" "<<tsos.localPosition().y()<<" "<<tsos.localPosition().z();
845  */
846  if(tsos.isValid()){
847  //check the angle of this tsos
848  float mom_x=tsos.localDirection().x();
849  float mom_y=tsos.localDirection().y();
850  float mom_z=tsos.localDirection().z();
851  //we took LOCAL momentum, i.e. respect to surface. Thus the plane is z=0
852  float angle=TMath::ASin(TMath::Abs(mom_z) / sqrt(pow(mom_x,2)+pow(mom_y,2)+pow(mom_z,2) ) );
853  if(!rejectLowAngleHits_ || angle>=TrackAngleCut_) angle_ok=true;// keep this hit
854  // else std::cout<<"Hit rejected because angle is "<< angle<<" ( <"<<TrackAngleCut_<<" )"<<std::endl;
855 
856  if(angle_ok && PXLcorrClusChargeCut_>0.0){
857  //
858  //get the hit from the TM and check that it is in the pixel
860  if(hitpointer->isValid()){
861  const TrackingRecHit *hit=(*hitpointer).hit();
862  if(GeomDetEnumerators::isTrackerPixel(theGeometry->geomDetSubDetector(hit->geographicalId().subdetId()))) {//do it only for pixel hits
863  corrcharge_ok=false;
864  float clust_alpha= atan2( mom_z, mom_x );
865  float clust_beta= atan2( mom_z, mom_y );
866 
867  //Now get the cluster charge
868 
869  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
870  float clust_charge=pixelhit->cluster()->charge();
871  float corr_clust_charge=clust_charge * sqrt( 1.0 / ( 1.0/pow( tan(clust_alpha), 2 ) +
872  1.0/pow( tan(clust_beta ), 2 ) +
873  1.0 )
874  );
875  //std::cout<<"xxxxx "<<clust_charge<<" "<<corr_clust_charge<<" " <<pixelhit->qBin()<<" "<<pixelhit->clusterProbability(1)<<" "<<pixelhit->clusterProbability(2)<< std::endl;
876  if(corr_clust_charge>PXLcorrClusChargeCut_){
877  corrcharge_ok=true;
878  }
879  } //end if hit is in pixel
880  }//end if hit is valid
881 
882  }//check corr cluster charge for pixel hits
883 
884  }//end if TSOS is valid
885  else{
886  edm::LogWarning("TrackerTrackHitFilter") <<"TSOS not valid ! Impossible to calculate track angle.";
887  }
888 
889  return angle_ok&&corrcharge_ok;
890 }//end TrackerTrackHitFilter::checkHitAngle
ConstRecHitPointer const & recHit() const
LocalVector localDirection() const
T y() const
Definition: PV3DBase.h:63
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
T Abs(T a)
Definition: MathUtil.h:49
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
edm::ESHandle< TrackerGeometry > theGeometry
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
TrajectoryStateOnSurface const & updatedState() const
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Our base class.
Definition: SiPixelRecHit.h:23
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
bool reco::modules::TrackerTrackHitFilter::checkPXLCorrClustCharge ( const TrajectoryMeasurement meas)
private

Definition at line 893 of file TrackerTrackHitFilter.cc.

References SiPixelRecHit::cluster(), TrackingRecHit::geographicalId(), GeomDetEnumerators::isTrackerStrip(), TrajectoryStateOnSurface::isValid(), TrajectoryStateOnSurface::localDirection(), funct::pow(), PXLcorrClusChargeCut_, TrajectoryMeasurement::recHit(), mathSSE::sqrt(), DetId::subdetId(), funct::tan(), theGeometry, TrajectoryMeasurement::updatedState(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

893  {
894  /*
895  Code taken from DPGAnalysis/SiPixelTools/plugins/PixelNtuplizer_RealData.cc
896  */
897 
898  bool corrcharge_ok=false;
899  //get the hit from the TM and check that it is in the pixel
901  if(!hitpointer->isValid()) return corrcharge_ok;
902  const TrackingRecHit *hit=(*hitpointer).hit();
903  if(GeomDetEnumerators::isTrackerStrip(theGeometry->geomDetSubDetector(hit->geographicalId().subdetId()))) {//SiStrip hit, skip
904  return corrcharge_ok;
905  }
906 
908  if(tsos.isValid()){
909  float mom_x=tsos.localDirection().x();
910  float mom_y=tsos.localDirection().y();
911  float mom_z=tsos.localDirection().z();
912  float clust_alpha= atan2( mom_z, mom_x );
913  float clust_beta= atan2( mom_z, mom_y );
914 
915  //Now get the cluster charge
916 
917  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
918  float clust_charge=pixelhit->cluster()->charge();
919  float corr_clust_charge=clust_charge * sqrt( 1.0 / ( 1.0/pow( tan(clust_alpha), 2 ) +
920  1.0/pow( tan(clust_beta ), 2 ) +
921  1.0 )
922  );
923  if(corr_clust_charge>PXLcorrClusChargeCut_)corrcharge_ok=true;
924 
925  }//end if TSOS is valid
926  return corrcharge_ok;
927 
928 }//end TrackerTrackHitFilter::checkPXLCorrClustCharge
ConstRecHitPointer const & recHit() const
LocalVector localDirection() const
T y() const
Definition: PV3DBase.h:63
bool isTrackerStrip(const GeomDetEnumerators::SubDetector m)
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
edm::ESHandle< TrackerGeometry > theGeometry
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
TrajectoryStateOnSurface const & updatedState() const
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Our base class.
Definition: SiPixelRecHit.h:23
bool reco::modules::TrackerTrackHitFilter::checkStoN ( const edm::EventSetup iSetup,
const DetId id,
const TrackingRecHit therechit 
)
private

Definition at line 739 of file TrackerTrackHitFilter.cc.

References checkPXLQuality_, SiStripRecHit1D::cluster(), SiStripRecHit2D::cluster(), SiPixelRecHit::clusterProbability(), Exception, SiPixelRecHit::hasFilledProb(), GeomDetEnumerators::isTrackerPixel(), GeomDetEnumerators::isTrackerStrip(), pxlTPLProbXY_, pxlTPLProbXYQ_, pxlTPLqBin_, SiPixelRecHit::qBin(), SiStripClusterInfo::signalOverNoise(), subdetStoN_, subdetStoNhighcut_, subdetStoNlowcut_, and theGeometry.

Referenced by checkHit().

739  {
740 
741  bool keepthishit=true;
742  // const uint32_t& recHitDetId = id.rawId();
743 
744  //check StoN only if subdet was set by the user
745  // int subdet_cnt=0;
746  int subdet_cnt=id.subdetId();
747 
748  // for(subdet_cnt=0;subdet_cnt<6; ++subdet_cnt){
749 
750  // if( subdetStoN_[subdet_cnt-1]&& (id.subdetId()==subdet_cnt) ){//check that hit is in a det belonging to a subdet where we decided to apply a S/N cut
751 
752 
753  if(GeomDetEnumerators::isTrackerStrip(theGeometry->geomDetSubDetector(id.subdetId()))) {
754  if( subdetStoN_[subdet_cnt-1]){//check that hit is in a det belonging to a subdet where we decided to apply a S/N cut
755  const std::type_info &type = typeid(*therechit);
756  const SiStripCluster* cluster;
757  if (type == typeid(SiStripRecHit2D)) {
758  const SiStripRecHit2D* hit = dynamic_cast<const SiStripRecHit2D*>(therechit);
759  if (hit!=0) cluster = &*(hit->cluster());
760  else{
761  edm::LogError("TrackerTrackHitFilter")<< "TrackerTrackHitFilter::checkStoN : Unknown valid tracker hit in subdet " << id.subdetId()<< "(detID="<<id.rawId()<<")\n ";
762  keepthishit = false;
763  }
764  }
765  else if (type == typeid(SiStripRecHit1D)) {
766  const SiStripRecHit1D* hit = dynamic_cast<const SiStripRecHit1D*>(therechit);
767  if (hit!=0) cluster = &*(hit->cluster());
768  else{
769  edm::LogError("TrackerTrackHitFilter")<< "TrackerTrackHitFilter::checkStoN : Unknown valid tracker hit in subdet " << id.subdetId()<< "(detID="<<id.rawId()<<")\n ";
770  keepthishit = false;
771  }
772  }
773  //the following two cases should not happen anymore since CMSSW > 2_0_X because of hit splitting in stereo modules
774  //const SiStripMatchedRecHit2D* matchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(therechit);
775  //const ProjectedSiStripRecHit2D* unmatchedhit = dynamic_cast<const ProjectedSiStripRecHit2D*>(therechit);
776  else{
777  throw cms::Exception("Unknown RecHit Type") << "RecHit of type " << type.name() <<
778  " not supported. (use c++filt to demangle the name)";
779  }
780 
781  if(keepthishit){
782  SiStripClusterInfo clusterInfo = SiStripClusterInfo( *cluster, iSetup, id.rawId());
783  if ( (subdetStoNlowcut_[subdet_cnt-1]>0) && (clusterInfo.signalOverNoise() < subdetStoNlowcut_[subdet_cnt-1]) ) keepthishit = false;
784  if ( (subdetStoNhighcut_[subdet_cnt-1]>0) && (clusterInfo.signalOverNoise() > subdetStoNhighcut_[subdet_cnt-1]) ) keepthishit = false;
785  //if(!keepthishit)std::cout<<"Hit rejected because of bad S/N: "<<clusterInfo.signalOverNoise()<<std::endl;
786  }
787 
788  }//end if subdetStoN_[subdet_cnt]&&...
789 
790  }//end if subdet_cnt >2
791  else if (GeomDetEnumerators::isTrackerPixel(theGeometry->geomDetSubDetector(id.subdetId()))){//pixel
792  //pixels have naturally a very low noise (because of their low capacitance). So the S/N cut is
793  //irrelevant in this case. Leave it dummy
794  keepthishit = true;
795 
796  /**************
797  * Cut on cluster charge corr by angle embedded in the checkHitAngle() function
798  *
799  *************/
800 
801  if(checkPXLQuality_){
802  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(therechit);
803  if(pixelhit!=0){
804  //std::cout << "ClusterCharge=" <<std::flush<<pixelhit->cluster()->charge() << std::flush;
805  float xyprob=pixelhit->clusterProbability(0);//x-y combined log_e probability of the pixel cluster
806  //singl x- and y-prob not stored sicne CMSSW 3_9_0
807  float xychargeprob =pixelhit->clusterProbability(1);//xy-prob * charge prob
808  // float chargeprob =pixelhit->clusterProbability(2);//charge prob
809  bool haspassed_tplreco= pixelhit->hasFilledProb(); //the cluster was associted to a template
810  int qbin =pixelhit->qBin(); //R==meas_charge/pred_charge: Qbin=0 ->R>1.5 , =1->1<R<1.5 ,=2->0.85<R<1 ,
811  // Qbin=3->0.95*Qminpred<R<0.85 ,=4->, =5->meas_charge<0.95*Qminpred
812 
813  // if(haspassed_tplreco) std::cout<<" CLUSTPROB=\t"<<xprob<<"\t"<<yprob<<"\t"<<combprob<<"\t"<<qbin<<std::endl;
814  // else std::cout<<"CLUSTPROBNOTDEF=\t"<<xprob<<"\t"<<yprob<<"\t"<<combprob<<"\t"<<qbin<<std::endl;
815 
816  keepthishit = false;
817  // std::cout<<"yyyyy "<<qbin<<" "<<xprob<<" "<<yprob<<std::endl;
818  if( haspassed_tplreco && xyprob>pxlTPLProbXY_ && xychargeprob>pxlTPLProbXYQ_ && qbin>pxlTPLqBin_[0] && qbin<=pxlTPLqBin_[1] )keepthishit = true;
819 
820  }
821  else {edm::LogInfo("TrackerTrackHitFilter")<<"HIT IN PIXEL ("<<subdet_cnt <<") but PixelRecHit is EMPTY!!!";}
822  }//end if check pixel quality flag
823  }
824  // else throw cms::Exception("TrackerTrackHitFilter") <<"Loop over subdetector out of range when applying the S/N cut: "<<subdet_cnt;
825 
826  // }//end loop on subdets
827 
828 
829  return keepthishit;
830 }//end CheckStoN
ClusterRef cluster() const
type
Definition: HCALResponse.h:21
bool hasFilledProb() const
float clusterProbability(unsigned int flags=0) const
bool isTrackerStrip(const GeomDetEnumerators::SubDetector m)
float signalOverNoise() const
ClusterRef cluster() const
int qBin() const
edm::ESHandle< TrackerGeometry > theGeometry
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
Our base class.
Definition: SiPixelRecHit.h:23
int reco::modules::TrackerTrackHitFilter::layerFromId ( const DetId id,
const TrackerTopology tTopo 
) const
private

Definition at line 932 of file TrackerTrackHitFilter.cc.

References TrackerTopology::layer().

Referenced by produceFromTrajectory().

933 {
934  return tTopo->layer(id);
935 }
unsigned int layer(const DetId &id) const
TrackCandidate reco::modules::TrackerTrackHitFilter::makeCandidate ( const reco::Track tk,
std::vector< TrackingRecHit * >::iterator  hitsBegin,
std::vector< TrackingRecHit * >::iterator  hitsEnd 
)
private

Definition at line 523 of file TrackerTrackHitFilter.cc.

References alongMomentum, anyDirection, Exception, reco::Track::innerDetId(), reco::Track::innerPosition(), trajectoryStateTransform::innerStateOnSurface(), reco::TrackBase::momentum(), reco::Track::outerDetId(), reco::Track::outerPosition(), trajectoryStateTransform::outerStateOnSurface(), trajectoryStateTransform::persistentState(), edm::OwnVector< T, P >::push_back(), edm::OwnVector< T, P >::reserve(), fileCollector::seed, reco::Track::seedDirection(), reco::Track::seedRef(), theGeometry, and theMagField.

Referenced by produce().

523  {
524 
525 
527  PTrajectoryStateOnDet state;
528  if ( pdir == anyDirection ) throw cms::Exception("UnimplementedFeature") << "Cannot work with tracks that have 'anyDirecton' \n";
529 
530  // double innerP=sqrt( pow(tk.innerMomentum().X(),2)+pow(tk.innerMomentum().Y(),2)+pow(tk.innerMomentum().Z(),2) );
531  // if ( (pdir == alongMomentum) == ( innerP >= tk.outerP() ) ) {
532 
533  if ( (pdir == alongMomentum) == ( (tk.outerPosition()-tk.innerPosition()).Dot(tk.momentum()) >= 0 ) ) {
534  // use inner state
536  state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(tk.innerDetId()) );
537  } else {
538  // use outer state
540  state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(tk.outerDetId()) );
541  }
544  ownHits.reserve(hitsEnd - hitsBegin);
545  for ( ; hitsBegin != hitsEnd; ++hitsBegin) {
546  //if(! (*hitsBegin)->isValid() ) std::cout<<"Putting in the trackcandidate an INVALID HIT !"<<std::endl;
547  ownHits.push_back( *hitsBegin );
548  }
549 
550  TrackCandidate cand(ownHits, seed, state, tk.seedRef());
551 
552  return cand;
553 }
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
PropagationDirection
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:662
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
void push_back(D *&d)
Definition: OwnVector.h:280
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:94
Definition: DetId.h:18
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:213
edm::ESHandle< TrackerGeometry > theGeometry
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:204
edm::ESHandle< MagneticField > theMagField
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
void reserve(size_t)
Definition: OwnVector.h:274
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
void reco::modules::TrackerTrackHitFilter::parseStoN ( const std::string &  str)
private

Definition at line 206 of file TrackerTrackHitFilter.cc.

References Exception, plotBeamSpotDB::first, match(), edm::second(), AlCaHLTBitMon_QueryRunRegistry::string, subdetStoN_, subdetStoNhighcut_, and subdetStoNlowcut_.

Referenced by TrackerTrackHitFilter().

206  {
207  // match a set of capital case chars (preceded by an arbitrary number of leading blanks),
208  //followed b an arbitrary number of blanks, one or more digits (not necessary, they cannot also be,
209  // another set of blank spaces and, again another *eventual* digit
210  // static boost::regex rule("\\s+([A-Z]+)(\\s+(\\d+)(\\.)?(\\d+))?(\\s+(\\d+)(\\.)?(\\d+))?");
211  static const boost::regex rule("([A-Z]+)"
212  "\\s*(\\d+\\.*\\d*)?"
213  "\\s*(\\d+\\.*\\d*)?");
214 
215 
216  boost::cmatch match;
217  std::string match_1;
218  std::string match_2;
219  std::string match_3;
220  // match and check it works
221  if (!regex_match(str.c_str(), match, rule)) {
222  throw cms::Exception("Configuration") << "Rule for S to N cut '" << str << "' not understood.\n";
223  }
224  else{
225  std::string match_0=match[0].second;
226  match_1=match[1].second;
227  match_2=match[2].second;
228  match_3=match[3].second;
229 
230  }
231 
232  int cnt=0;
233  float subdet_ind[6];
234  for(cnt=0;cnt<6;cnt++ ){
235  subdet_ind[cnt]=-1.0;
236  }
237 
238 
239  bool doALL=false;
240  std::string match_1a(match[1].first,match[1].second);
241  if (strncmp(match[1].first, "ALL", 3) == 0) doALL=true;
242  if (doALL || strncmp(match[1].first, "PXB", 3) == 0) subdet_ind[0] = +1.0;
243  if (doALL || strncmp(match[1].first, "PXE", 3) == 0) subdet_ind[1] = +1.0;
244  if (doALL || strncmp(match[1].first, "TIB", 3) == 0) subdet_ind[2] = +1.0;
245  if (doALL || strncmp(match[1].first, "TID", 3) == 0) subdet_ind[3] = +1.0;
246  if (doALL || strncmp(match[1].first, "TOB", 3) == 0) subdet_ind[4] = +1.0;
247  if (doALL || strncmp(match[1].first, "TEC", 3) == 0) subdet_ind[5] = +1.0;
248 
249  for(cnt=0;cnt<6;cnt++ ){//loop on subdets
250  if(subdet_ind[cnt]>0.0){
251  subdetStoN_[cnt]=true;
252  if (match[2].first != match[2].second) {
253  subdetStoNlowcut_[cnt] = atof(match[2].first);
254  }
255  if (match[3].first != match[3].second ) {
256  subdetStoNhighcut_[cnt] = atof(match[3].first);
257  }
258  edm::LogInfo("TrackerTrackHitFilter") <<"Setting thresholds*&^ for subdet #"<<cnt+1<<" = "<<subdetStoNlowcut_[cnt]<<" - "<<subdetStoNhighcut_[cnt];
259  }
260  }
261 
262  bool correct_regex=false;
263  for(cnt=0;cnt<6;cnt++ ){//check that the regex was correct
264  if(subdetStoN_[cnt])correct_regex=true;
265  }
266 
267  if (!correct_regex) {
268  throw cms::Exception("Configuration") << "Detector '" << match_1a << "' not understood in parseStoN. Should be PXB, PXE, TIB, TID, TOB, TEC.\n";
269  }
270 
271  //std::cout<<"Reached end of parseStoN"<<std::endl;
272 }//end parseStoN
U second(std::pair< T, U > const &p)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
void reco::modules::TrackerTrackHitFilter::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Implements edm::stream::EDProducerBase.

Definition at line 359 of file TrackerTrackHitFilter.cc.

References begin, edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > >::const_iterator, end, edm::EventID::event(), edm::false, edm::EventSetup::get(), edm::Event::getByToken(), edm::EventBase::id(), iEvt, iRun, edm::Ref< C, T, F >::key(), makeCandidate(), minimumHits_, convertSQLitetoXML_cfg::output, produceFromTrack(), produceFromTrajectory(), edm::Event::put(), replaceWithInactiveHits_, edm::EventID::run(), stripBackInvalidHits_, stripFrontInvalidHits_, theGeometry, theMagField, tokenTracks, tokenTrajTrack, testEve_cfg::tracks, and useTrajectories_.

360 {
361  //Dump Run and Event
362  iRun=iEvent.id().run();
363  iEvt=iEvent.id().event();
364 
365  // read with View, so we can read also a TrackRefVector
368 
369  if(useTrajectories_)iEvent.getByToken(tokenTrajTrack, assoMap);
370  else iEvent.getByToken(tokenTracks, tracks);
371 
372  // read from EventSetup
375  //iSetup.get<TransientRecHitRecord>().get("WithTrackAngle",theBuilder);
376  // RHBuilder= theBuilder.product();
377 
378 
379  // prepare output collection
380  size_t candcollsize;
381  if(useTrajectories_)candcollsize= assoMap->size();
382  else candcollsize=tracks->size();
383  std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection());
384 
385  output->reserve(candcollsize);
386 
387  // working area and tools
388  std::vector<TrackingRecHit *> hits;
389 
390  if(useTrajectories_){
391  for (TrajTrackAssociationCollection::const_iterator itass = assoMap->begin(); itass != assoMap->end(); ++itass){
392 
393  const edm::Ref<std::vector<Trajectory> >traj = itass->key;//trajectory in the collection
394  const reco::TrackRef tkref = itass->val;//associated track track in the collection
395  //std::cout<<"The hit collection has size "<<hits.size()<<" (should be 0) while the track contains initially "<< tkref->recHitsEnd() - tkref->recHitsBegin()<<std::endl;
396 
397  const Track *trk = &(*tkref);
398  const Trajectory * myTrajectory= &(*traj);
399  produceFromTrajectory(iSetup,myTrajectory,hits);
400 
401  std::vector<TrackingRecHit *>::iterator begin = hits.begin(), end = hits.end();
402 
403  // strip invalid hits at the beginning
405  while ( (begin != end) && ( (*begin)->isValid() == false ) ) ++begin;
406  }
407  // strip invalid hits at the end
408  if (stripBackInvalidHits_ && (begin != end)) {
409  --end;
410  while ( (begin != end) && ( (*end)->isValid() == false ) ) --end;
411  ++end;
412  }
413 
414  // if we still have some hits build the track candidate
416  int nvalidhits=0;
417  for( std::vector<TrackingRecHit *>::iterator ithit=begin;ithit!=end;++ithit){
418  if( (*ithit)->isValid())nvalidhits++;
419  }
420  if(nvalidhits >= int(minimumHits_)){
421  output->push_back( makeCandidate ( *trk, begin, end ) );
422  }
423  }
424  else{//all invalid hits have been already kicked out
425  if ((end - begin) >= int(minimumHits_)) {
426  output->push_back( makeCandidate ( *trk, begin, end ) );
427  }
428  }
429 
430 
431  // if we still have some hits build the track candidate
432  //if ((end - begin) >= int(minimumHits_)) {
433  // output->push_back( makeCandidate ( *trk, begin, end ) );
434  //}
435 
436 
437  // now delete the hits not used by the candidate
438  for (begin = hits.begin(), end = hits.end(); begin != end; ++begin) {
439  if (*begin) delete *begin;
440  }
441  hits.clear();
442  } // loop on trajectories
443 
444 
445  }
446  else{ //use plain tracks
447 
448  // loop on tracks
449  for (std::vector<reco::Track>::const_iterator ittrk = tracks->begin(), edtrk = tracks->end(); ittrk != edtrk; ++ittrk) {
450 
451  // std::cout<<"The hit collection has size "<<hits.size()<<" (should be 0) while the track contains initially "<< ittrk->recHitsEnd() - ittrk->recHitsBegin()<<std::endl;
452 
453  const Track *trk = &(*ittrk);
454 
455  produceFromTrack(iSetup,trk,hits);
456  //-----------------------
457  /*
458  std::cout<<"Hit collection in output has size "<<hits.size()<<". Dumping hit positions..."<<std::endl;
459  for (std::vector<TrackingRecHit *>::iterator ith = hits.begin(), edh = hits.end(); ith != edh; ++ith) {
460  const TrackingRecHit *myhit = *(ith);
461  TransientTrackingRecHit::RecHitPointer ttrh;
462  float radius=0.0;
463  float xx=-999.0,yy=-999.0,zz=-999.0;
464  unsigned int myid=0;
465  if(myhit->isValid()){
466  ttrh = RHBuilder->build(myhit);
467  xx=ttrh->globalPosition().x();
468  yy=ttrh->globalPosition().y();
469  zz=ttrh->globalPosition().z();
470  radius = sqrt( pow(xx,2)+pow(yy,2) );
471  myid=myhit->geographicalId().rawId();
472  }
473  std::cout<<"-$-$ OUTPUT Hit position: ( "<<xx<<" , " <<yy<<" , " <<zz<<" ) , RADIUS = " <<radius<<" on DetID= "<< myid<<std::endl;
474  }//end loop on hits
475  */
476  //-----------------------
477 
478 
479  std::vector<TrackingRecHit *>::iterator begin = hits.begin(), end = hits.end();
480  // std::cout << "Back in the main producer (TRK), the final hit collection has size " << hits.size() << std::endl;
481  // strip invalid hits at the beginning
483  while ( (begin != end) && ( (*begin)->isValid() == false ) ) ++begin;
484  }
485  // strip invalid hits at the end
486  if (stripBackInvalidHits_ && (begin != end)) {
487  --end;
488  while ( (begin != end) && ( (*end)->isValid() == false ) ) --end;
489  ++end;
490  }
491 
492  // if we still have some hits build the track candidate
494  int nvalidhits=0;
495  for( std::vector<TrackingRecHit *>::iterator ithit=begin;ithit!=end;++ithit){
496  if( (*ithit)->isValid())nvalidhits++;
497  }
498  if(nvalidhits >= int(minimumHits_)){
499  output->push_back( makeCandidate ( *ittrk, begin, end ) );
500  }
501 
502  }
503  else{//all invalid hits have been already kicked out
504  if ((end - begin) >= int(minimumHits_)) {
505  output->push_back( makeCandidate ( *ittrk, begin, end ) );
506  }
507  }
508 
509  // now delete the hits not used by the candidate
510  for (begin = hits.begin(), end = hits.end(); begin != end; ++begin) {
511  if (*begin) delete *begin;
512  }
513  hits.clear();
514  } // loop on tracks
515  }//end else useTracks
516 
517  // std::cout<<"OUTPUT SIZE: "<<output->size()<<std::endl;
518 
519  iEvent.put(output);
520 }
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< TrackCandidate > TrackCandidateCollection
edm::EDGetTokenT< reco::TrackCollection > tokenTracks
key_type key() const
Accessor for product key.
Definition: Ref.h:264
void produceFromTrajectory(const edm::EventSetup &iSetup, const Trajectory *itt, std::vector< TrackingRecHit * > &hits)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
void produceFromTrack(const edm::EventSetup &iSetup, const Track *itt, std::vector< TrackingRecHit * > &hits)
#define end
Definition: vmac.h:37
tuple tracks
Definition: testEve_cfg.py:39
edm::ESHandle< TrackerGeometry > theGeometry
const T & get() const
Definition: EventSetup.h:56
edm::EventID id() const
Definition: EventBase.h:60
#define begin
Definition: vmac.h:30
edm::EDGetTokenT< TrajTrackAssociationCollection > tokenTrajTrack
volatile std::atomic< bool > shutdown_flag false
TrackCandidate makeCandidate(const reco::Track &tk, std::vector< TrackingRecHit * >::iterator hitsBegin, std::vector< TrackingRecHit * >::iterator hitsEnd)
edm::ESHandle< MagneticField > theMagField
void reco::modules::TrackerTrackHitFilter::produceFromTrack ( const edm::EventSetup iSetup,
const Track itt,
std::vector< TrackingRecHit * > &  hits 
)

Definition at line 555 of file TrackerTrackHitFilter.cc.

References checkHit(), TrackingRecHit::clone(), TrackingRecHit::det(), cond::rpcobgas::detid, TrackingRecHit::geographicalId(), TrackingRecHit::inactive, TrackingRecHit::isValid(), DetId::rawId(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), replaceWithInactiveHits_, and stripAllInvalidHits_.

Referenced by produce().

555  {
556 
557  // loop on tracks
558  hits.clear(); // extra safety
559 
560  for (trackingRecHit_iterator ith = itt->recHitsBegin(), edh = itt->recHitsEnd(); ith != edh; ++ith) {
561  const TrackingRecHit * hit = (*ith); // ith is an iterator on edm::Ref to rechit
562 
563  DetId detid = hit->geographicalId();
564 
565  //check that the hit is a real hit and not a constraint
566  if(hit->isValid() && hit==0 && detid.rawId()==0) continue;
567 
568  int verdict=checkHit(iSetup,detid,hit);
569  if (verdict == 0) {
570  // just copy the hit
571  hits.push_back(hit->clone());
572  }
573  else if(verdict<-2){//hit rejected because did not pass the selections
574  // still, if replaceWithInactiveHits is true we have to put a new hit
576  hits.push_back(new InvalidTrackingRecHit(*(hit->det()), TrackingRecHit::inactive));
577  }
578  }
579  else if(verdict==-2) hits.push_back(hit->clone());//hit not in the tracker
580  else if(verdict==-1){ //hit not valid
581  if (!stripAllInvalidHits_) {
582  hits.push_back(hit->clone());
583  }
584  }
585  } // loop on hits
586 
587 }//end TrackerTrackHitFilter::produceFromTrack
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const GeomDet * det() const
virtual TrackingRecHit * clone() const =0
Definition: DetId.h:18
int checkHit(const edm::EventSetup &iSetup, const DetId &detid, const TrackingRecHit *hit)
bool isValid() const
DetId geographicalId() const
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
void reco::modules::TrackerTrackHitFilter::produceFromTrajectory ( const edm::EventSetup iSetup,
const Trajectory itt,
std::vector< TrackingRecHit * > &  hits 
)

—OverlapEnd

—OverlapBegin

—OverlapEnd

Definition at line 590 of file TrackerTrackHitFilter.cc.

References checkHit(), checkHitAngle(), TrackingRecHit::clone(), TrackingRecHit::det(), cond::rpcobgas::detid, TrackingRecHit::geographicalId(), edm::EventSetup::get(), TrackingRecHit::inactive, GeomDetEnumerators::isTrackerPixel(), GeomDetEnumerators::isTrackerStrip(), layerFromId(), Trajectory::measurements(), nOverlaps, edm::ESHandle< class >::product(), rejectLowAngleHits_, replaceWithInactiveHits_, sideFromId(), SiStripDetId::stereo(), stripAllInvalidHits_, DetId::subdetId(), tagOverlaps_, and theGeometry.

Referenced by produce().

590  {
591  hits.clear(); // extra safety
592  nOverlaps=0;
593 
594  //Retrieve tracker topology from geometry
596  iSetup.get<TrackerTopologyRcd>().get(tTopoHand);
597  const TrackerTopology *tTopo=tTopoHand.product();
598 
599 
600  std::vector<TrajectoryMeasurement> tmColl =itt->measurements();
601 
602  //---OverlapBegin needed eventually for overlaps, but I must create them here in any case
603  const TrajectoryMeasurement* previousTM(0);
604  DetId previousId(0);
605  //int previousLayer(-1);
607 
608  int constrhits=0;
609 
610  for(std::vector<TrajectoryMeasurement>::const_iterator itTrajMeas = tmColl.begin(); itTrajMeas!=tmColl.end(); itTrajMeas++){
611  TransientTrackingRecHit::ConstRecHitPointer hitpointer = itTrajMeas->recHit();
612 
613  //check that the hit is a real hit and not a constraint
614  if(hitpointer->isValid() && hitpointer->hit()==0){constrhits++; continue;}
615 
616  const TrackingRecHit *hit=((*hitpointer).hit());
617  DetId detid = hit->geographicalId();
618  int verdict=checkHit(iSetup,detid,hit);
619 
620  if (verdict == 0) {
621  if( rejectLowAngleHits_ && !checkHitAngle(*itTrajMeas) ){//check angle of track on module if requested
622  verdict=-6;//override previous verdicts
623  }
624  }
625 
626  /*
627  //this has been included in checkHitAngle(*itTrajMeas)
628  if (verdict == 0) {
629  if( PXLcorrClusChargeCut_>0.0 && !checkPXLCorrClustCharge(*itTrajMeas) ){//check angle of track on module if requested
630  verdict=-7;//override previous verdicts
631  }
632  }
633  */
634 
635 
636  if(verdict==0){// Hit TAKEN !!!!
637 
638  if(tagOverlaps_){
639  //std::cout<<"Looking for overlaps in Run="<<iRun<<" , Event ="<<iEvt<<std::flush;
640 
641  int side(sideFromId(detid,tTopo));//side 0=barrel, 1=minus , 2=plus
642  int layer(layerFromId(detid,tTopo));//layer or disk
643  int subDet = detid.subdetId();
644  //std::cout << " Check Subdet #" <<subDet << ", layer = " <<layer<<" stereo: "<< ((subDet > 2)?(SiStripDetId(detid).stereo()):2);
645 
646  if ( ( previousTM!=0 )&& (layer!=-1 )) {
647  //std::cout<<"A previous TM exists! "<<std::endl;
648  for (std::vector<TrajectoryMeasurement>::const_iterator itmCompare =itTrajMeas-1;itmCompare >= tmColl.begin() && itmCompare > itTrajMeas - 4;--itmCompare){
649 
650  DetId compareId = itmCompare->recHit()->geographicalId();
651  if ( subDet != compareId.subdetId() ||
652  side != sideFromId(compareId,tTopo) ||
653  layer != layerFromId(compareId,tTopo)) break;
654  if (!itmCompare->recHit()->isValid()) continue;
655  if(GeomDetEnumerators::isTrackerPixel(theGeometry->geomDetSubDetector(detid.subdetId())) ||
656  (GeomDetEnumerators::isTrackerStrip(theGeometry->geomDetSubDetector(detid.subdetId())) &&
657  SiStripDetId(detid).stereo()==SiStripDetId(compareId).stereo()))
658  {//if either pixel or strip stereo module
659  // overlapHits.push_back(std::make_pair(&(*itmCompare),&(*itm)));
660  //std::cout<< "Adding pair "<< ((subDet >2)?(SiStripDetId(detid).stereo()):2)
661  // << " from SubDet = "<<subDet<<" , layer = " << layer<<" Run:"<<iRun<<"\tEv: "<<iEvt<<"\tId1: "<<compareId.rawId()<<"\tId2: "<<detid.rawId()<<std::endl;
662  // if(abs(compareId.rawId()-detid.rawId())==1)std::cout<<"These two are from the same det! Id1= "<<detid.rawId()<<" has stereo type "<<SiStripDetId(detid).stereo() <<"\tId2: "<<compareId.rawId()<<" has stereo type "<<SiStripDetId(compareId).stereo()<<std::endl;
664  // if(detid.rawId()<compareId.rawId()){
665  // std::cout<< "+++ "<< "\t"<<iRun<<"\t"<<iEvt<<"\t"<<detid.rawId()<<"\t"<<compareId.rawId()<<std::endl;
666  // }
667  //else std::cout<< "+++ "<< "\t"<<iRun<<"\t"<<iEvt<<"\t"<<compareId.rawId()<<"\t"<<detid.rawId()<<std::endl;
668 
669  nOverlaps++;
670  break;
671  }
672  }//end second loop on TM for overlap tagging
673 
674  }//end if ( (layer!=-1 )&&(acceptLayer[subDet]))
675 
676  previousTM = &(* itTrajMeas);
677  previousId = detid;
678  //previousLayer = layer;
679  }//end if look for overlaps
681 
682  hits.push_back(hit->clone()); //just copy it
683  }//end if HIT TAKEN
684  else if(verdict<-2){//hit rejected because did not pass the selections
685  // still, if replaceWithInactiveHits is true we have to put a new hit
687  hits.push_back(new InvalidTrackingRecHit(*hit->det(), TrackingRecHit::inactive));
688  }
689  }
690  else if(verdict==-2) hits.push_back(hit->clone());//hit not in the tracker
691  else if(verdict==-1){ //hit not valid
692  if (!stripAllInvalidHits_) {
693  hits.push_back(hit->clone());
694  }
695  }
696  } // loop on hits
697 
698  std::reverse(hits.begin(),hits.end());
699  // std::cout<<"Finished producefromTrajecotries. Nhits in final coll"<<hits.size() <<" Nconstraints="<<constrhits<<std::endl;
700 } //end TrackerTrackHitFilter::produceFromTrajectories
uint32_t stereo() const
Definition: SiStripDetId.h:176
bool isTrackerStrip(const GeomDetEnumerators::SubDetector m)
DataContainer const & measurements() const
Definition: Trajectory.h:203
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
const GeomDet * det() const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
virtual TrackingRecHit * clone() const =0
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
Definition: DetId.h:18
int layerFromId(const DetId &id, const TrackerTopology *tTopo) const
int checkHit(const edm::EventSetup &iSetup, const DetId &detid, const TrackingRecHit *hit)
edm::ESHandle< TrackerGeometry > theGeometry
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
int sideFromId(const DetId &id, const TrackerTopology *tTopo) const
bool checkHitAngle(const TrajectoryMeasurement &meas)
DetId geographicalId() const
int reco::modules::TrackerTrackHitFilter::sideFromId ( const DetId id,
const TrackerTopology tTopo 
) const
private

Definition at line 937 of file TrackerTrackHitFilter.cc.

References TrackerTopology::side().

Referenced by produceFromTrajectory().

938 {
939  return tTopo->side(id);
940 }
unsigned int side(const DetId &id) const

Member Data Documentation

bool reco::modules::TrackerTrackHitFilter::checkPXLQuality_
private

Definition at line 137 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN().

std::string reco::modules::TrackerTrackHitFilter::CMNSubtractionMode_
private

Definition at line 123 of file TrackerTrackHitFilter.cc.

std::vector<uint32_t> reco::modules::TrackerTrackHitFilter::detsToIgnore_
private

Definition at line 130 of file TrackerTrackHitFilter.cc.

Referenced by checkHit(), and TrackerTrackHitFilter().

edm::EventNumber_t reco::modules::TrackerTrackHitFilter::iEvt
private

Definition at line 113 of file TrackerTrackHitFilter.cc.

Referenced by produce().

edm::RunNumber_t reco::modules::TrackerTrackHitFilter::iRun
private

Definition at line 112 of file TrackerTrackHitFilter.cc.

Referenced by produce().

size_t reco::modules::TrackerTrackHitFilter::minimumHits_
private

Definition at line 115 of file TrackerTrackHitFilter.cc.

Referenced by produce().

int reco::modules::TrackerTrackHitFilter::nOverlaps
private

Definition at line 154 of file TrackerTrackHitFilter.cc.

Referenced by produceFromTrajectory().

double reco::modules::TrackerTrackHitFilter::PXLcorrClusChargeCut_
private
double reco::modules::TrackerTrackHitFilter::pxlTPLProbXY_
private

Definition at line 138 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN().

double reco::modules::TrackerTrackHitFilter::pxlTPLProbXYQ_
private

Definition at line 139 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN().

std::vector<int32_t> reco::modules::TrackerTrackHitFilter::pxlTPLqBin_
private

Definition at line 140 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN(), and TrackerTrackHitFilter().

bool reco::modules::TrackerTrackHitFilter::rejectBadStoNHits_
private

Definition at line 122 of file TrackerTrackHitFilter.cc.

Referenced by checkHit(), and TrackerTrackHitFilter().

bool reco::modules::TrackerTrackHitFilter::rejectLowAngleHits_
private
bool reco::modules::TrackerTrackHitFilter::replaceWithInactiveHits_
private
std::vector<Rule> reco::modules::TrackerTrackHitFilter::rules_
private

Definition at line 131 of file TrackerTrackHitFilter.cc.

Referenced by checkHit(), and TrackerTrackHitFilter().

edm::InputTag reco::modules::TrackerTrackHitFilter::src_
private

Definition at line 110 of file TrackerTrackHitFilter.cc.

Referenced by TrackerTrackHitFilter().

bool reco::modules::TrackerTrackHitFilter::stripAllInvalidHits_
private
bool reco::modules::TrackerTrackHitFilter::stripBackInvalidHits_
private

Definition at line 119 of file TrackerTrackHitFilter.cc.

Referenced by produce().

bool reco::modules::TrackerTrackHitFilter::stripFrontInvalidHits_
private

Definition at line 118 of file TrackerTrackHitFilter.cc.

Referenced by produce().

std::vector<bool> reco::modules::TrackerTrackHitFilter::subdetStoN_
private

Definition at line 124 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN(), parseStoN(), and TrackerTrackHitFilter().

std::vector<double> reco::modules::TrackerTrackHitFilter::subdetStoNhighcut_
private

Definition at line 126 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN(), parseStoN(), and TrackerTrackHitFilter().

std::vector<double> reco::modules::TrackerTrackHitFilter::subdetStoNlowcut_
private

Definition at line 125 of file TrackerTrackHitFilter.cc.

Referenced by checkStoN(), parseStoN(), and TrackerTrackHitFilter().

bool reco::modules::TrackerTrackHitFilter::tagOverlaps_
private

Definition at line 153 of file TrackerTrackHitFilter.cc.

Referenced by produceFromTrajectory().

edm::ESHandle<TransientTrackingRecHitBuilder> reco::modules::TrackerTrackHitFilter::theBuilder
private

Definition at line 148 of file TrackerTrackHitFilter.cc.

edm::ESHandle<TrackerGeometry> reco::modules::TrackerTrackHitFilter::theGeometry
private
edm::ESHandle<MagneticField> reco::modules::TrackerTrackHitFilter::theMagField
private

Definition at line 146 of file TrackerTrackHitFilter.cc.

Referenced by makeCandidate(), and produce().

edm::EDGetTokenT<reco::TrackCollection> reco::modules::TrackerTrackHitFilter::tokenTracks
private

Definition at line 150 of file TrackerTrackHitFilter.cc.

Referenced by produce(), and TrackerTrackHitFilter().

edm::EDGetTokenT<TrajTrackAssociationCollection> reco::modules::TrackerTrackHitFilter::tokenTrajTrack
private

Definition at line 151 of file TrackerTrackHitFilter.cc.

Referenced by produce(), and TrackerTrackHitFilter().

double reco::modules::TrackerTrackHitFilter::TrackAngleCut_
private

Definition at line 135 of file TrackerTrackHitFilter.cc.

Referenced by checkHitAngle(), and TrackerTrackHitFilter().

bool reco::modules::TrackerTrackHitFilter::useTrajectories_
private

Definition at line 133 of file TrackerTrackHitFilter.cc.

Referenced by produce(), and TrackerTrackHitFilter().