#include <FastElectronSeedGenerator.h>
Class to generate the trajectory seed from two hits in the pixel detector which have been found compatible with an ECAL cluster.
Description: Top algorithm producing ElectronSeeds, ported from FAMOS
Implementation: future redesign...
Definition at line 51 of file FastElectronSeedGenerator.h.
Definition at line 57 of file FastElectronSeedGenerator.h.
Definition at line 56 of file FastElectronSeedGenerator.h.
Definition at line 59 of file FastElectronSeedGenerator.h.
Definition at line 58 of file FastElectronSeedGenerator.h.
Definition at line 62 of file FastElectronSeedGenerator.h.
FastElectronSeedGenerator::FastElectronSeedGenerator | ( | const edm::ParameterSet & | pset, |
double | pTMin, | ||
const edm::InputTag & | beamSpot | ||
) |
Definition at line 53 of file FastElectronSeedGenerator.cc.
References gather_cfg::cout, fromTrackerSeeds_, edm::ParameterSet::getParameter(), myGSPixelMatcher, and searchInTIDTEC.
: dynamicphiroad_(pset.getParameter<bool>("dynamicPhiRoad")), lowPtThreshold_(pset.getParameter<double>("LowPtThreshold")), highPtThreshold_(pset.getParameter<double>("HighPtThreshold")), sizeWindowENeg_(pset.getParameter<double>("SizeWindowENeg")), phimin2_(pset.getParameter<double>("PhiMin2")), phimax2_(pset.getParameter<double>("PhiMax2")), deltaPhi1Low_(pset.getParameter<double>("DeltaPhi1Low")), deltaPhi1High_(pset.getParameter<double>("DeltaPhi1High")), deltaPhi2_(pset.getParameter<double>("DeltaPhi2")), pTMin2(pTMin*pTMin), myGSPixelMatcher(0), fromTrackerSeeds_(pset.getParameter<bool>("fromTrackerSeeds")), theUpdator(0), thePropagator(0), // theMeasurementTracker(0), // theNavigationSchool(0) theSetup(0), theBeamSpot(beamSpot), pts_(0) { #ifdef FAMOS_DEBUG std::cout << "FromTrackerSeeds = " << fromTrackerSeeds_ << std::endl; #endif // Instantiate the pixel hit matcher searchInTIDTEC = pset.getParameter<bool>("searchInTIDTEC"); myGSPixelMatcher = new FastPixelHitMatcher(pset.getParameter<double>("ePhiMin1"), pset.getParameter<double>("ePhiMax1"), pset.getParameter<double>("pPhiMin1"), pset.getParameter<double>("pPhiMax1"), pset.getParameter<double>("PhiMin2"), pset.getParameter<double>("PhiMax2"), pset.getParameter<double>("z2MinB"), pset.getParameter<double>("z2MaxB"), pset.getParameter<double>("r2MinF"), pset.getParameter<double>("r2MaxF"), pset.getParameter<double>("rMinI"), pset.getParameter<double>("rMaxI"), pset.getParameter<bool>("searchInTIDTEC")); }
FastElectronSeedGenerator::~FastElectronSeedGenerator | ( | ) |
Definition at line 98 of file FastElectronSeedGenerator.cc.
References myGSPixelMatcher, thePropagator, and theUpdator.
{ // delete theNavigationSchool; delete myGSPixelMatcher; // delete myMatchPos; delete thePropagator; delete theUpdator; }
void FastElectronSeedGenerator::addASeedToThisCluster | ( | edm::Ref< reco::SuperClusterCollection > | seedCluster, |
std::vector< TrackerRecHit > & | theHits, | ||
const TrajectorySeed & | theTrackerSeed, | ||
std::vector< reco::ElectronSeed > & | result | ||
) | [private] |
Definition at line 308 of file FastElectronSeedGenerator.cc.
References alongMomentum, BSPosition_, FastPixelHitMatcher::compatibleHits(), gather_cfg::cout, deltaPhi1High_, deltaPhi1Low_, deltaPhi2_, dir, dynamicphiroad_, fromTrackerSeeds_, FastPixelHitMatcher::getVertex(), highPtThreshold_, lowPtThreshold_, myGSPixelMatcher, prepareElTrackSeed(), pts_, recHits_, asciidump::s, FastPixelHitMatcher::set1stLayer(), FastPixelHitMatcher::set2ndLayer(), reco::ElectronSeed::setCaloCluster(), funct::sin(), sizeWindowENeg_, v, and TrackValidation_HighPurity_cff::valid.
Referenced by run().
{ float clusterEnergy = seedCluster->energy(); GlobalPoint clusterPos(seedCluster->position().x(), seedCluster->position().y(), seedCluster->position().z()); const GlobalPoint vertexPos(BSPosition_.x(),BSPosition_.y(),BSPosition_.z()); #ifdef FAMOS_DEBUG std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] " << "new supercluster with energy: " << clusterEnergy << std::endl; std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] " << "and position: " << clusterPos << std::endl; std::cout << "Vertex position : " << vertexPos << std::endl; #endif //Here change the deltaPhi window of the first pixel layer in function of the seed pT if (dynamicphiroad_) { float clusterEnergyT = clusterEnergy*sin(seedCluster->position().theta()) ; float deltaPhi1 = 0.875/clusterEnergyT + 0.055; if (clusterEnergyT < lowPtThreshold_) deltaPhi1= deltaPhi1Low_; if (clusterEnergyT > highPtThreshold_) deltaPhi1= deltaPhi1High_; float ephimin1 = -deltaPhi1*sizeWindowENeg_ ; float ephimax1 = deltaPhi1*(1.-sizeWindowENeg_); float pphimin1 = -deltaPhi1*(1.-sizeWindowENeg_); float pphimax1 = deltaPhi1*sizeWindowENeg_; float phimin2 = -deltaPhi2_/2.; float phimax2 = deltaPhi2_/2.; myGSPixelMatcher->set1stLayer(ephimin1,ephimax1,pphimin1,pphimax1); myGSPixelMatcher->set2ndLayer(phimin2,phimax2); } PropagationDirection dir = alongMomentum; // Find the best pixel pair compatible with the cluster std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> > compatPixelHits = myGSPixelMatcher->compatibleHits(clusterPos, vertexPos, clusterEnergy, theHits); // The corresponding origin vertex double vertexZ = myGSPixelMatcher->getVertex(); GlobalPoint theVertex(BSPosition_.x(),BSPosition_.y(),vertexZ); // Create the Electron pixel seed. if (!compatPixelHits.empty() ) { #ifdef FAMOS_DEBUG std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] " << " electron compatible hits found " << std::endl; #endif // Pixel-matching case: create the seed from scratch if (!fromTrackerSeeds_) { std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> >::iterator v; for (v = compatPixelHits.begin(); v != compatPixelHits.end(); ++v ) { bool valid = prepareElTrackSeed(v->first,v->second, theVertex); if (valid) { reco::ElectronSeed s(*pts_,recHits_,dir) ; s.setCaloCluster(reco::ElectronSeed::CaloClusterRef(seedCluster)) ; result.push_back(s); delete pts_; pts_=0; } } // Here we take instead the seed from a-priori seed collection } else { reco::ElectronSeed s(theTrackerSeed); s.setCaloCluster(reco::ElectronSeed::CaloClusterRef(seedCluster)) ; result.push_back(s); } } #ifdef FAMOS_DEBUG else std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] " << " No electron compatible hits found " << std::endl; #endif // And return ! return ; }
bool FastElectronSeedGenerator::prepareElTrackSeed | ( | ConstRecHitPointer | outerhit, |
ConstRecHitPointer | innerhit, | ||
const GlobalPoint & | vertexPos | ||
) | [private] |
Definition at line 406 of file FastElectronSeedGenerator.cc.
References edm::OwnVector< T, P >::clear(), TrajectoryStateOnSurface::isValid(), LogDebug, TrajectoryStateTransform::persistentState(), PropagatorWithMaterial::propagate(), pts_, edm::OwnVector< T, P >::push_back(), recHits_, FreeTrajectoryState::setCurvilinearError(), thePropagator, theSetup, and transformer_.
Referenced by addASeedToThisCluster().
{ // debug prints LogDebug("") <<"[FastElectronSeedGenerator::prepareElTrackSeed] " << "inner PixelHit x,y,z "<<innerhit->globalPosition(); LogDebug("") <<"[FastElectronSeedGenerator::prepareElTrackSeed] " << "outer PixelHit x,y,z "<<outerhit->globalPosition(); pts_=0; // make a spiral from the two hits and the vertex position FastHelix helix(outerhit->globalPosition(),innerhit->globalPosition(),vertexPos,*theSetup); if ( !helix.isValid()) return false; FreeTrajectoryState fts = helix.stateAtVertex(); // Give infinite errors to start the fit (no pattern recognition here). AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID(); fts.setCurvilinearError(errorMatrix*100.); TrajectoryStateOnSurface propagatedState = thePropagator->propagate(fts,innerhit->det()->surface()) ; if (!propagatedState.isValid()) return false; // The persitent trajectory state pts_ = transformer_.persistentState(propagatedState, innerhit->geographicalId().rawId()); // The corresponding rechits recHits_.clear(); recHits_.push_back(innerhit->hit()->clone()); recHits_.push_back(outerhit->hit()->clone()); return true; }
void FastElectronSeedGenerator::run | ( | edm::Event & | e, |
const reco::SuperClusterRefVector & | sclRefs, | ||
const SiTrackerGSMatchedRecHit2DCollection * | theGSRecHits, | ||
const edm::SimTrackContainer * | theSimTracks, | ||
TrajectorySeedCollection * | seeds, | ||
reco::ElectronSeedCollection & | out | ||
) |
Definition at line 142 of file FastElectronSeedGenerator.cc.
References addASeedToThisCluster(), BSPosition_, gather_cfg::cout, fromTrackerSeeds_, edm::RangeMap< ID, C, P >::get(), edm::Event::getByLabel(), i, edm::EventBase::id(), edm::RangeMap< ID, C, P >::ids(), TrackerRecHit::layerNumber(), LogDebug, CoreSimTrack::momentum(), myGSPixelMatcher, pTMin2, TrajectorySeed::recHits(), TrackerRecHit::ringNumber(), searchInTIDTEC, FastPixelHitMatcher::set1stLayerZRange(), findQualityFiles::size, edm::RefVector< C, T, F >::size(), edm::RangeMap< ID, C, P >::size(), mathSSE::sqrt(), TrackerRecHit::subDetId(), theBeamSpot, theInitialSeedColl, and theTrackerGeometry.
Referenced by FastElectronSeedProducer::produce().
{ // Take the seed collection. theInitialSeedColl=seeds; // Get the beam spot edm::Handle<reco::BeamSpot> recoBeamSpotHandle; e.getByLabel(theBeamSpot,recoBeamSpotHandle); // Get its position BSPosition_ = recoBeamSpotHandle->position(); double sigmaZ=recoBeamSpotHandle->sigmaZ(); double sigmaZ0Error=recoBeamSpotHandle->sigmaZ0Error(); double sq=sqrt(sigmaZ*sigmaZ+sigmaZ0Error*sigmaZ0Error); double zmin1 = BSPosition_.z()-3*sq; double zmax1 = BSPosition_.z()+3*sq; #ifdef FAMOS_DEBUG std::cout << "Z Range for pixel matcher : " << zmin1 << " " << BSPosition_.z() << " " << zmax1 << std::endl; #endif myGSPixelMatcher->set1stLayerZRange(zmin1,zmax1); // A map of vector of pixel seeds, for each clusters std::map<unsigned,std::vector<reco::ElectronSeed> > myPixelSeeds; // No seeding attempted if no hits ! if(theGSRecHits->size() == 0) return; if ( !fromTrackerSeeds_ ) { // The vector of simTrack Id's carrying GSRecHits const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids(); // Loop over the simTrack carrying GSRecHits for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) { unsigned simTrackId = theSimTrackIds[tkId]; const SimTrack& theSimTrack = (*theSimTracks)[simTrackId]; // Request a minimum pT for the sim track if ( theSimTrack.momentum().perp2() < pTMin2 ) continue; // Request a minimum number of RecHits (total and in the pixel detector) unsigned numberOfRecHits = 0; // The vector of rechits for seeding // 1) Cluster-pixel match seeding: // Save a collection of Pixel +TEC +TID hits for seeding electrons std::vector<unsigned> layerHit(6,static_cast<unsigned>(0)); // const SiTrackerGSMatchedRecHit2D *hit; TrackerRecHit currentHit; std::vector<TrackerRecHit> theHits; TrajectorySeed theTrackerSeed; SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId); SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first; SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second; SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit; SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit2; for ( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) { ++numberOfRecHits; currentHit = TrackerRecHit(&(*iterRecHit),theTrackerGeometry); if ( ( currentHit.subDetId() <= 2 ) || // Pixel Hits // Add TID/TEC (optional) ( searchInTIDTEC && ( ( currentHit.subDetId() == 3 && currentHit.ringNumber() < 3 && currentHit.layerNumber() < 3 ) || // TID first two rings, first two layers ( currentHit.subDetId() == 6 && currentHit.ringNumber() < 3 && currentHit.layerNumber() < 3 ) ) ) ) // TEC first two rings, first two layers theHits.push_back(currentHit); } // At least 3 hits if ( numberOfRecHits < 3 ) continue; // At least 2 pixel hits if ( theHits.size() < 2 ) continue; // Loop over clusters unsigned csize = sclRefs.size(); for (unsigned int i=0;i<csize;++i) { // Find the pixel seeds (actually only the best one is returned) LogDebug ("run") << "new cluster, calling addAseedFromThisCluster"; addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]); } } // 2) Check if the seed is in the a-priori seed collection } else { // Loop over the tracker seed #ifdef FAMOS_DEBUG std::cout << "We have " << seeds->size() << " tracker seeds!" << std::endl; #endif for (unsigned int i=0;i<seeds->size();++i) { TrackerRecHit currentHit; std::vector<TrackerRecHit> theHits; const TrajectorySeed& theTrackerSeed = (*seeds)[i]; TrajectorySeed::range theSeedRange=theTrackerSeed.recHits(); TrajectorySeed::const_iterator theSeedRangeIteratorBegin = theSeedRange.first; TrajectorySeed::const_iterator theSeedRangeIteratorEnd = theSeedRange.second; TrajectorySeed::const_iterator theSeedItr = theSeedRangeIteratorBegin; for ( ; theSeedItr != theSeedRangeIteratorEnd; ++theSeedItr ) { const SiTrackerGSMatchedRecHit2D * theSeedingRecHit = (const SiTrackerGSMatchedRecHit2D*) (&(*theSeedItr)); currentHit = TrackerRecHit(theSeedingRecHit,theTrackerGeometry); theHits.push_back(currentHit); } // Loop over clusters unsigned csize = sclRefs.size(); for (unsigned int i=0;i<csize;++i) { // Find the pixel seeds (actually only the best one is returned) #ifdef FAMOS_DEBUG std::cout << "new cluster, calling addAseedFromThisCluster" << std::endl; #endif addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]); } // End loop over clusters } // End loop over seeds } // end else // Back to the expected collection std::map<unsigned,std::vector<reco::ElectronSeed> >::const_iterator px = myPixelSeeds.begin(); std::map<unsigned,std::vector<reco::ElectronSeed> >::const_iterator pxEnd = myPixelSeeds.end(); for ( ; px!=pxEnd; ++px ) { unsigned nSeeds = (px->second).size(); for ( unsigned ipx = 0; ipx<nSeeds; ++ipx ) { out.push_back((px->second)[ipx]); reco::ElectronSeed is = px->second[ipx]; } } LogDebug ("run") << ": For event "<<e.id(); LogDebug ("run") <<"Nr of superclusters: "<<sclRefs.size() <<", no. of ElectronSeeds found = " << out.size(); #ifdef FAMOS_DEBUG std::cout << ": For event "<<e.id() << std::endl; std::cout <<"Nr of superclusters: "<<sclRefs.size() <<", no. of ElectronSeeds found = " << out.size() << std::endl; #endif }
void FastElectronSeedGenerator::setup | ( | bool | ) |
Referenced by setupES().
void FastElectronSeedGenerator::setupES | ( | const edm::EventSetup & | setup | ) |
Definition at line 109 of file FastElectronSeedGenerator.cc.
References alongMomentum, geometry, edm::EventSetup::get(), myGSPixelMatcher, FastPixelHitMatcher::setES(), setup(), theGeomSearchTracker, theMagField, theMagneticFieldMap, thePropagator, theSetup, theTrackerGeometry, and theTrackerInteractionGeometry.
Referenced by FastElectronSeedProducer::beginRun().
{ theSetup= &setup; edm::ESHandle<MagneticField> pMF; setup.get<IdealMagneticFieldRecord>().get(pMF); theMagField = &(*pMF); edm::ESHandle<TrackerGeometry> geometry; setup.get<TrackerDigiGeometryRecord>().get(geometry); theTrackerGeometry = &(*geometry); edm::ESHandle<GeometricSearchTracker> recoGeom; setup.get<TrackerRecoGeometryRecord>().get( recoGeom ); theGeomSearchTracker = &(*recoGeom); edm::ESHandle<TrackerInteractionGeometry> interGeom; setup.get<TrackerInteractionGeometryRecord>().get( interGeom ); theTrackerInteractionGeometry = &(*interGeom); edm::ESHandle<MagneticFieldMap> fieldMap; setup.get<MagneticFieldMapRecord>().get(fieldMap); theMagneticFieldMap = &(*fieldMap); thePropagator = new PropagatorWithMaterial(alongMomentum,.000511,&(*theMagField)); myGSPixelMatcher->setES(theMagneticFieldMap, theTrackerGeometry, theGeomSearchTracker, theTrackerInteractionGeometry); }
Definition at line 111 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster(), and run().
float FastElectronSeedGenerator::deltaPhi1High_ [private] |
Definition at line 105 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
float FastElectronSeedGenerator::deltaPhi1Low_ [private] |
Definition at line 105 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
float FastElectronSeedGenerator::deltaPhi2_ [private] |
Definition at line 106 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
bool FastElectronSeedGenerator::dynamicphiroad_ [private] |
Definition at line 99 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
bool FastElectronSeedGenerator::fromTrackerSeeds_ [private] |
Definition at line 117 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster(), FastElectronSeedGenerator(), and run().
float FastElectronSeedGenerator::highPtThreshold_ [private] |
Definition at line 102 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
float FastElectronSeedGenerator::lowPtThreshold_ [private] |
Definition at line 101 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
Definition at line 113 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster(), FastElectronSeedGenerator(), run(), setupES(), and ~FastElectronSeedGenerator().
float FastElectronSeedGenerator::phimax2_ [private] |
Definition at line 104 of file FastElectronSeedGenerator.h.
float FastElectronSeedGenerator::phimin2_ [private] |
Definition at line 104 of file FastElectronSeedGenerator.h.
double FastElectronSeedGenerator::pTMin2 [private] |
Definition at line 110 of file FastElectronSeedGenerator.h.
Referenced by run().
Definition at line 132 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster(), and prepareElTrackSeed().
Definition at line 131 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster(), and prepareElTrackSeed().
bool FastElectronSeedGenerator::searchInTIDTEC [private] |
Definition at line 107 of file FastElectronSeedGenerator.h.
Referenced by FastElectronSeedGenerator(), and run().
float FastElectronSeedGenerator::sizeWindowENeg_ [private] |
Definition at line 103 of file FastElectronSeedGenerator.h.
Referenced by addASeedToThisCluster().
const edm::InputTag FastElectronSeedGenerator::theBeamSpot [private] |
Definition at line 129 of file FastElectronSeedGenerator.h.
Referenced by run().
const GeometricSearchTracker* FastElectronSeedGenerator::theGeomSearchTracker [private] |
Definition at line 122 of file FastElectronSeedGenerator.h.
Referenced by setupES().
Definition at line 116 of file FastElectronSeedGenerator.h.
Referenced by run().
const MagneticField* FastElectronSeedGenerator::theMagField [private] |
Definition at line 119 of file FastElectronSeedGenerator.h.
Referenced by setupES().
const MagneticFieldMap* FastElectronSeedGenerator::theMagneticFieldMap [private] |
Definition at line 120 of file FastElectronSeedGenerator.h.
Referenced by setupES().
Definition at line 126 of file FastElectronSeedGenerator.h.
Referenced by prepareElTrackSeed(), setupES(), and ~FastElectronSeedGenerator().
const edm::EventSetup* FastElectronSeedGenerator::theSetup [private] |
Definition at line 128 of file FastElectronSeedGenerator.h.
Referenced by prepareElTrackSeed(), and setupES().
const TrackerGeometry* FastElectronSeedGenerator::theTrackerGeometry [private] |
Definition at line 121 of file FastElectronSeedGenerator.h.
const TrackerInteractionGeometry* FastElectronSeedGenerator::theTrackerInteractionGeometry [private] |
Definition at line 123 of file FastElectronSeedGenerator.h.
Referenced by setupES().
KFUpdator* FastElectronSeedGenerator::theUpdator [private] |
Definition at line 125 of file FastElectronSeedGenerator.h.
Referenced by ~FastElectronSeedGenerator().
Definition at line 130 of file FastElectronSeedGenerator.h.
Referenced by prepareElTrackSeed().
double FastElectronSeedGenerator::zmax1_ [private] |
Definition at line 109 of file FastElectronSeedGenerator.h.
double FastElectronSeedGenerator::zmin1_ [private] |
Definition at line 109 of file FastElectronSeedGenerator.h.