CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTracker/SpecialSeedGenerators/src/BeamHaloPairGenerator.cc

Go to the documentation of this file.
00001 #include "RecoTracker/SpecialSeedGenerators/interface/BeamHaloPairGenerator.h"
00002 #include "RecoTracker/TkSeedingLayers/interface/SeedingLayerSetsBuilder.h"
00003 typedef TransientTrackingRecHit::ConstRecHitPointer SeedingHit;
00004 
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 using namespace ctfseeding;
00007 
00008 
00009 BeamHaloPairGenerator::BeamHaloPairGenerator(const edm::ParameterSet& conf): conf_(conf){
00010         edm::LogInfo("CtfSpecialSeedGenerator|BeamHaloPairGenerator") << "Constructing BeamHaloPairGenerator";
00011         theMaxTheta=conf.getParameter<double>("maxTheta");
00012         theMaxTheta=fabs(sin(theMaxTheta));
00013 } 
00014 
00015 
00016 SeedingLayerSets BeamHaloPairGenerator::init(const edm::EventSetup& es){
00017         edm::ParameterSet leyerPSet = conf_.getParameter<edm::ParameterSet>("LayerPSet");
00018         SeedingLayerSetsBuilder lsBuilder(leyerPSet);
00019         SeedingLayerSets lss = lsBuilder.layers(es);
00020         return lss;     
00021 }
00022 
00023 
00024 const OrderedSeedingHits& BeamHaloPairGenerator::run(const TrackingRegion& region,
00025                                                     const edm::Event& e,
00026                                                     const edm::EventSetup& es){
00027         hitPairs.clear();
00028         SeedingLayerSets lss = init(es);
00029         SeedingLayerSets::const_iterator iLss;
00030         for (iLss = lss.begin(); iLss != lss.end(); iLss++){
00031                 SeedingLayers ls = *iLss;
00032                 if (ls.size() != 2){
00033                         throw cms::Exception("CtfSpecialSeedGenerator") << "You are using " << ls.size() <<" layers in set instead of 2 ";
00034                 }       
00035                 std::vector<SeedingHit> innerHits  = region.hits(e, es, &ls[0]);
00036                 std::vector<SeedingHit> outerHits  = region.hits(e, es, &ls[1]);
00037                 std::vector<SeedingHit>::const_iterator iOuterHit;
00038                 for (iOuterHit = outerHits.begin(); iOuterHit != outerHits.end(); iOuterHit++){
00039                         std::vector<SeedingHit>::const_iterator iInnerHit;
00040                         for (iInnerHit = innerHits.begin(); iInnerHit != innerHits.end(); iInnerHit++){
00041                           //do something in there... if necessary
00042                           const TransientTrackingRecHit::ConstRecHitPointer & crhpi = *iInnerHit;
00043                           const TransientTrackingRecHit::ConstRecHitPointer & crhpo =  *iOuterHit;
00044                           GlobalVector d=crhpo->globalPosition() - crhpi->globalPosition();
00045                           double ABSsinDtheta = fabs(sin(d.theta()));
00046                           LogDebug("BeamHaloPairGenerator")<<"position1: "<<crhpo->globalPosition()
00047                                                            <<" position2: "<<crhpi->globalPosition()
00048                                                            <<" |sin(Dtheta)|: "<< ABSsinDtheta <<((ABSsinDtheta>theMaxTheta)?" skip":" keep");
00049 
00050                           if (ABSsinDtheta>theMaxTheta) {;continue;}
00051 
00052                           hitPairs.push_back(OrderedHitPair(*iInnerHit,
00053                                                             *iOuterHit));
00054                         }
00055                 }
00056         }
00057         return hitPairs;
00058 }