CMS 3D CMS Logo

ElectronGSPixelSeedGenerator Class Reference

Class to generate the trajectory seed from two hits in the pixel detector which have been found compatible with an ECAL cluster. More...

#include <FastSimulation/EgammaElectronAlgos/interface/ElectronGSPixelSeedGenerator.h>

List of all members.

Public Types

typedef
TransientTrackingRecHit::ConstRecHitPointer 
ConstRecHitPointer
enum  mode { HLT, offline, unknown }
typedef edm::OwnVector
< TrackingRecHit
PRecHitContainer
typedef
TransientTrackingRecHit::RecHitContainer 
RecHitContainer
typedef
TransientTrackingRecHit::RecHitPointer 
RecHitPointer

Public Member Functions

 ElectronGSPixelSeedGenerator (const edm::ParameterSet &pset, double pTMin, const edm::InputTag &beamSpot)
void run (edm::Event &e, const reco::SuperClusterRefVector &sclRefs, const SiTrackerGSMatchedRecHit2DCollection *theGSRecHits, const edm::SimTrackContainer *theSimTracks, TrajectorySeedCollection *seeds, reco::ElectronPixelSeedCollection &out)
void setup (bool)
void setupES (const edm::EventSetup &setup)
 ~ElectronGSPixelSeedGenerator ()

Private Member Functions

void addASeedToThisCluster (edm::Ref< reco::SuperClusterCollection > seedCluster, std::vector< TrackerRecHit > &theHits, const TrajectorySeed &theTrackerSeed, std::vector< reco::ElectronPixelSeed > &result)
bool prepareElTrackSeed (ConstRecHitPointer outerhit, ConstRecHitPointer innerhit, const GlobalPoint &vertexPos)

Private Attributes

math::XYZPoint BSPosition_
float deltaPhi1High_
float deltaPhi1Low_
float deltaPhi2_
bool dynamicphiroad_
bool fromTrackerSeeds_
float highPtThreshold_
float lowPtThreshold_
GSPixelHitMatchermyGSPixelMatcher
float phimax2_
float phimin2_
double pTMin2
PTrajectoryStateOnDetpts_
PRecHitContainer recHits_
bool searchInTIDTEC
float sizeWindowENeg_
const edm::InputTag theBeamSpot
const GeometricSearchTrackertheGeomSearchTracker
TrajectorySeedCollectiontheInitialSeedColl
const MagneticFieldtheMagField
const MagneticFieldMaptheMagneticFieldMap
PropagatorWithMaterialthePropagator
const edm::EventSetuptheSetup
const TrackerGeometrytheTrackerGeometry
const TrackerInteractionGeometrytheTrackerInteractionGeometry
KFUpdatortheUpdator
TrajectoryStateTransform transformer_
double zmax1_
double zmin1_


Detailed Description

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 ElectronPixelSeeds, ported from FAMOS.

Author:
Patrick Janot
Implementation: future redesign...

Definition at line 51 of file ElectronGSPixelSeedGenerator.h.


Member Typedef Documentation

typedef TransientTrackingRecHit::ConstRecHitPointer ElectronGSPixelSeedGenerator::ConstRecHitPointer

Definition at line 57 of file ElectronGSPixelSeedGenerator.h.

typedef edm::OwnVector<TrackingRecHit> ElectronGSPixelSeedGenerator::PRecHitContainer

Definition at line 56 of file ElectronGSPixelSeedGenerator.h.

typedef TransientTrackingRecHit::RecHitContainer ElectronGSPixelSeedGenerator::RecHitContainer

Definition at line 59 of file ElectronGSPixelSeedGenerator.h.

typedef TransientTrackingRecHit::RecHitPointer ElectronGSPixelSeedGenerator::RecHitPointer

Definition at line 58 of file ElectronGSPixelSeedGenerator.h.


Member Enumeration Documentation

enum ElectronGSPixelSeedGenerator::mode

Enumerator:
HLT 
offline 
unknown 

Definition at line 62 of file ElectronGSPixelSeedGenerator.h.

00062 {HLT, offline, unknown};  //to be used later


Constructor & Destructor Documentation

ElectronGSPixelSeedGenerator::ElectronGSPixelSeedGenerator ( const edm::ParameterSet pset,
double  pTMin,
const edm::InputTag beamSpot 
)

Definition at line 53 of file ElectronGSPixelSeedGenerator.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), fromTrackerSeeds_, edm::ParameterSet::getParameter(), myGSPixelMatcher, and searchInTIDTEC.

00057   :   
00058   dynamicphiroad_(pset.getParameter<bool>("dynamicPhiRoad")),
00059   lowPtThreshold_(pset.getParameter<double>("LowPtThreshold")),
00060   highPtThreshold_(pset.getParameter<double>("HighPtThreshold")),
00061   sizeWindowENeg_(pset.getParameter<double>("SizeWindowENeg")),
00062   phimin2_(pset.getParameter<double>("PhiMin2")),      
00063   phimax2_(pset.getParameter<double>("PhiMax2")),
00064   deltaPhi1Low_(pset.getParameter<double>("DeltaPhi1Low")),
00065   deltaPhi1High_(pset.getParameter<double>("DeltaPhi1High")),
00066   deltaPhi2_(pset.getParameter<double>("DeltaPhi2")),
00067   pTMin2(pTMin*pTMin),
00068   myGSPixelMatcher(0),
00069   fromTrackerSeeds_(pset.getParameter<bool>("fromTrackerSeeds")),
00070   theUpdator(0), thePropagator(0), 
00071   //   theMeasurementTracker(0), 
00072   //   theNavigationSchool(0)
00073   theSetup(0), theBeamSpot(beamSpot), pts_(0)
00074 {
00075 
00076 #ifdef FAMOS_DEBUG
00077   std::cout << "FromTrackerSeeds  = " << fromTrackerSeeds_ << std::endl;
00078 #endif
00079 
00080   // Instantiate the pixel hit matcher
00081   searchInTIDTEC = pset.getParameter<bool>("searchInTIDTEC");
00082   myGSPixelMatcher = new GSPixelHitMatcher(pset.getParameter<double>("ePhiMin1"), 
00083                                            pset.getParameter<double>("ePhiMax1"),
00084                                            pset.getParameter<double>("pPhiMin1"),
00085                                            pset.getParameter<double>("pPhiMax1"),               
00086                                            pset.getParameter<double>("PhiMin2"),
00087                                            pset.getParameter<double>("PhiMax2"),
00088                                            pset.getParameter<double>("z2MinB"),
00089                                            pset.getParameter<double>("z2MaxB"),
00090                                            pset.getParameter<double>("r2MinF"),
00091                                            pset.getParameter<double>("r2MaxF"),
00092                                            pset.getParameter<double>("rMinI"),
00093                                            pset.getParameter<double>("rMaxI"),
00094                                            pset.getParameter<bool>("searchInTIDTEC"));
00095 
00096 }

ElectronGSPixelSeedGenerator::~ElectronGSPixelSeedGenerator (  ) 

Definition at line 98 of file ElectronGSPixelSeedGenerator.cc.

References myGSPixelMatcher, thePropagator, and theUpdator.

00098                                                             {
00099 
00100   //  delete theNavigationSchool;
00101   delete myGSPixelMatcher;
00102   //  delete myMatchPos;
00103   delete thePropagator;
00104   delete theUpdator;
00105 
00106 }


Member Function Documentation

void ElectronGSPixelSeedGenerator::addASeedToThisCluster ( edm::Ref< reco::SuperClusterCollection seedCluster,
std::vector< TrackerRecHit > &  theHits,
const TrajectorySeed theTrackerSeed,
std::vector< reco::ElectronPixelSeed > &  result 
) [private]

Definition at line 308 of file ElectronGSPixelSeedGenerator.cc.

References alongMomentum, BSPosition_, GSPixelHitMatcher::compatibleHits(), GenMuonPlsPt100GeV_cfg::cout, deltaPhi1High_, deltaPhi1Low_, deltaPhi2_, dir, dynamicphiroad_, lat::endl(), fromTrackerSeeds_, GSPixelHitMatcher::getVertex(), highPtThreshold_, lowPtThreshold_, myGSPixelMatcher, prepareElTrackSeed(), pts_, recHits_, s, GSPixelHitMatcher::set1stLayer(), GSPixelHitMatcher::set2ndLayer(), funct::sin(), sizeWindowENeg_, v, and TrackValidation_HighPurity_cff::valid.

Referenced by run().

00312 {
00313 
00314   float clusterEnergy = seedCluster->energy();
00315   GlobalPoint clusterPos(seedCluster->position().x(),
00316                          seedCluster->position().y(), 
00317                          seedCluster->position().z());
00318   const GlobalPoint vertexPos(BSPosition_.x(),BSPosition_.y(),BSPosition_.z());
00319 
00320 #ifdef FAMOS_DEBUG
00321   std::cout << "[ElectronGSPixelSeedGenerator::seedsFromThisCluster] " 
00322             << "new supercluster with energy: " << clusterEnergy << std::endl;
00323   std::cout << "[ElectronGSPixelSeedGenerator::seedsFromThisCluster] " 
00324             << "and position: " << clusterPos << std::endl;
00325   std::cout << "Vertex position : " << vertexPos << std::endl;
00326 #endif
00327 
00328   //Here change the deltaPhi window of the first pixel layer in function of the seed pT
00329   if (dynamicphiroad_) {
00330     float clusterEnergyT = clusterEnergy*sin(seedCluster->position().theta()) ;
00331     
00332     float deltaPhi1 = 0.875/clusterEnergyT + 0.055; 
00333     if (clusterEnergyT < lowPtThreshold_) deltaPhi1= deltaPhi1Low_;
00334     if (clusterEnergyT > highPtThreshold_) deltaPhi1= deltaPhi1High_;
00335     
00336     float ephimin1 = -deltaPhi1*sizeWindowENeg_ ;
00337     float ephimax1 =  deltaPhi1*(1.-sizeWindowENeg_);
00338     float pphimin1 = -deltaPhi1*(1.-sizeWindowENeg_);
00339     float pphimax1 =  deltaPhi1*sizeWindowENeg_;
00340     
00341     float phimin2  = -deltaPhi2_/2.;
00342     float phimax2  =  deltaPhi2_/2.;
00343     
00344     myGSPixelMatcher->set1stLayer(ephimin1,ephimax1,pphimin1,pphimax1);
00345     myGSPixelMatcher->set2ndLayer(phimin2,phimax2);
00346     
00347   }
00348 
00349 
00350 
00351   PropagationDirection dir = alongMomentum;
00352    
00353   // Find the best pixel pair compatible with the cluster
00354   std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> > compatPixelHits = 
00355     myGSPixelMatcher->compatibleHits(clusterPos, vertexPos, clusterEnergy, theHits);
00356     
00357   // The corresponding origin vertex
00358   double vertexZ = myGSPixelMatcher->getVertex();
00359   GlobalPoint theVertex(BSPosition_.x(),BSPosition_.y(),vertexZ); 
00360   
00361   // Create the Electron pixel seed.
00362   if (!compatPixelHits.empty() ) {
00363 #ifdef FAMOS_DEBUG
00364     std::cout << "[ElectronGSPixelSeedGenerator::seedsFromThisCluster] " 
00365               << " electron compatible hits found " << std::endl;
00366 #endif
00367     // Pixel-matching case: create the seed from scratch
00368     if (!fromTrackerSeeds_) {
00369       
00370       std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> >::iterator v;
00371       for (v = compatPixelHits.begin(); v != compatPixelHits.end(); ++v ) {
00372         
00373         bool valid = prepareElTrackSeed(v->first,v->second, theVertex);
00374         if (valid) {
00375           reco::ElectronPixelSeed s= reco::ElectronPixelSeed(seedCluster,*pts_,recHits_,dir);
00376           result.push_back(s);
00377           delete pts_;
00378           pts_=0;
00379         }
00380       }  
00381       
00382     // Here we take instead the seed from a-priori seed collection
00383     } else {
00384 
00385       reco::ElectronPixelSeed seed(seedCluster,theTrackerSeed);
00386       result.push_back(seed);
00387 
00388     }
00389 
00390   }
00391 #ifdef FAMOS_DEBUG
00392     else
00393       std::cout << "[ElectronGSPixelSeedGenerator::seedsFromThisCluster] " 
00394                 << " No electron compatible hits found " << std::endl;
00395 #endif
00396 
00397  
00398   // And return !
00399   return ;
00400 
00401 }

bool ElectronGSPixelSeedGenerator::prepareElTrackSeed ( ConstRecHitPointer  outerhit,
ConstRecHitPointer  innerhit,
const GlobalPoint vertexPos 
) [private]

Definition at line 403 of file ElectronGSPixelSeedGenerator.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().

00406 {
00407   
00408   // debug prints
00409   LogDebug("") <<"[ElectronGSPixelSeedGenerator::prepareElTrackSeed] " 
00410                << "inner PixelHit   x,y,z "<<innerhit->globalPosition();
00411   LogDebug("") <<"[ElectronGSPixelSeedGenerator::prepareElTrackSeed] " 
00412                << "outer PixelHit   x,y,z "<<outerhit->globalPosition();
00413 
00414   pts_=0;
00415   // make a spiral from the two hits and the vertex position
00416   FastHelix helix(outerhit->globalPosition(),innerhit->globalPosition(),vertexPos,*theSetup);
00417   if ( !helix.isValid()) return false;
00418 
00419   FreeTrajectoryState fts = helix.stateAtVertex();
00420 
00421   // Give infinite errors to start the fit (no pattern recognition here). 
00422   AlgebraicSymMatrix errorMatrix(5,1);
00423   fts.setCurvilinearError(errorMatrix*100.);
00424 
00425    TrajectoryStateOnSurface propagatedState = thePropagator->propagate(fts,innerhit->det()->surface()) ;
00426   if (!propagatedState.isValid()) return false;
00427 
00428   // The persitent trajectory state
00429   pts_ =  transformer_.persistentState(propagatedState, innerhit->geographicalId().rawId());
00430 
00431   // The corresponding rechits
00432   recHits_.clear();
00433   recHits_.push_back(innerhit->hit()->clone());
00434   recHits_.push_back(outerhit->hit()->clone());  
00435 
00436 
00437   return true;
00438 
00439 }

void ElectronGSPixelSeedGenerator::run ( edm::Event e,
const reco::SuperClusterRefVector sclRefs,
const SiTrackerGSMatchedRecHit2DCollection theGSRecHits,
const edm::SimTrackContainer theSimTracks,
TrajectorySeedCollection seeds,
reco::ElectronPixelSeedCollection out 
)

Definition at line 142 of file ElectronGSPixelSeedGenerator.cc.

References addASeedToThisCluster(), BSPosition_, GenMuonPlsPt100GeV_cfg::cout, lat::endl(), fromTrackerSeeds_, edm::RangeMap< ID, C, P >::get(), edm::Event::getByLabel(), i, edm::Event::id(), edm::RangeMap< ID, C, P >::ids(), TrackerRecHit::layerNumber(), LogDebug, CoreSimTrack::momentum(), myGSPixelMatcher, pTMin2, TrajectorySeed::recHits(), TrackerRecHit::ringNumber(), searchInTIDTEC, GSPixelHitMatcher::set1stLayerZRange(), simTrackId, edm::RefVector< C, T, F >::size(), edm::RangeMap< ID, C, P >::size(), size, funct::sqrt(), TrackerRecHit::subDetId(), theBeamSpot, theHits, theInitialSeedColl, and theTrackerGeometry.

Referenced by ElectronGSPixelSeedProducer::produce().

00147                                                                             {
00148 
00149   // Take the seed collection.
00150   theInitialSeedColl=seeds;
00151 
00152   // Get the beam spot
00153   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00154   e.getByLabel(theBeamSpot,recoBeamSpotHandle); 
00155 
00156   // Get its position
00157   BSPosition_ = recoBeamSpotHandle->position();
00158   double sigmaZ=recoBeamSpotHandle->sigmaZ();
00159   double sigmaZ0Error=recoBeamSpotHandle->sigmaZ0Error();
00160   double sq=sqrt(sigmaZ*sigmaZ+sigmaZ0Error*sigmaZ0Error);
00161   double zmin1 = BSPosition_.z()-3*sq;
00162   double zmax1 = BSPosition_.z()+3*sq;
00163 #ifdef FAMOS_DEBUG
00164   std::cout << "Z Range for pixel matcher : " << zmin1 << " " << BSPosition_.z() << " " << zmax1 << std::endl;
00165 #endif
00166   myGSPixelMatcher->set1stLayerZRange(zmin1,zmax1);
00167 
00168   // A map of vector of pixel seeds, for each clusters
00169   std::map<unsigned,std::vector<reco::ElectronPixelSeed> > myPixelSeeds;
00170 
00171   // No seeding attempted if no hits !
00172   if(theGSRecHits->size() == 0) return;    
00173   
00174   if ( !fromTrackerSeeds_ ) {
00175     
00176     // The vector of simTrack Id's carrying GSRecHits
00177     const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
00178     
00179     // Loop over the simTrack carrying GSRecHits
00180     for ( unsigned tkId=0;  tkId != theSimTrackIds.size(); ++tkId ) {
00181       
00182       unsigned simTrackId = theSimTrackIds[tkId];
00183       const SimTrack& theSimTrack = (*theSimTracks)[simTrackId]; 
00184       
00185       // Request a minimum pT for the sim track
00186       if ( theSimTrack.momentum().perp2() < pTMin2 ) continue;
00187       
00188       // Request a minimum number of RecHits (total and in the pixel detector)
00189       unsigned numberOfRecHits = 0;
00190       
00191       // The vector of rechits for seeding
00192       
00193       // 1) Cluster-pixel match seeding:
00194       //    Save a collection of Pixel +TEC +TID hits for seeding electrons
00195       std::vector<unsigned> layerHit(6,static_cast<unsigned>(0));
00196       // const SiTrackerGSMatchedRecHit2D *hit;
00197       TrackerRecHit currentHit;
00198       std::vector<TrackerRecHit> theHits;
00199       TrajectorySeed theTrackerSeed;
00200       
00201       SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
00202       SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00203       SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd   = theRecHitRange.second;
00204       SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00205       SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit2;
00206       
00207       for ( iterRecHit = theRecHitRangeIteratorBegin; 
00208             iterRecHit != theRecHitRangeIteratorEnd; 
00209             ++iterRecHit) { 
00210         ++numberOfRecHits;
00211         
00212         currentHit = TrackerRecHit(&(*iterRecHit),theTrackerGeometry);
00213         if ( ( currentHit.subDetId() <= 2 ) ||  // Pixel Hits
00214              // Add TID/TEC (optional)
00215              ( searchInTIDTEC && 
00216                ( ( currentHit.subDetId() == 3 && 
00217                    currentHit.ringNumber() < 3 && 
00218                    currentHit.layerNumber() < 3 ) || // TID first two rings, first two layers
00219                  ( currentHit.subDetId() == 6 && 
00220                    currentHit.ringNumber() < 3 && 
00221                    currentHit.layerNumber() < 3 ) ) ) ) // TEC first two rings, first two layers
00222           theHits.push_back(currentHit);
00223       }
00224 
00225       // At least 3 hits
00226       if ( numberOfRecHits < 3 ) continue;
00227       
00228       // At least 2 pixel hits
00229       if ( theHits.size() < 2 ) continue;
00230       
00231       // Loop over clusters 
00232       
00233       unsigned csize = sclRefs.size();
00234       for  (unsigned int i=0;i<csize;++i) {
00235         
00236         // Find the pixel seeds (actually only the best one is returned)
00237         LogDebug ("run") << "new cluster, calling addAseedFromThisCluster";
00238         addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]);
00239         
00240       }
00241     
00242     }
00243   // 2) Check if the seed is in the a-priori seed collection
00244   } else { 
00245 
00246     // Loop over the tracker seed
00247 #ifdef FAMOS_DEBUG
00248     std::cout << "We have " << seeds->size() << " tracker seeds!" << std::endl;
00249 #endif
00250     for (unsigned int i=0;i<seeds->size();++i) {
00251       
00252       TrackerRecHit currentHit;
00253       std::vector<TrackerRecHit> theHits;
00254       const TrajectorySeed& theTrackerSeed = (*seeds)[i];
00255       TrajectorySeed::range theSeedRange=theTrackerSeed.recHits();
00256       TrajectorySeed::const_iterator theSeedRangeIteratorBegin = theSeedRange.first;
00257       TrajectorySeed::const_iterator theSeedRangeIteratorEnd   = theSeedRange.second;
00258       TrajectorySeed::const_iterator theSeedItr = theSeedRangeIteratorBegin;
00259       
00260       for ( ; theSeedItr != theSeedRangeIteratorEnd; ++theSeedItr ) { 
00261         const SiTrackerGSMatchedRecHit2D * theSeedingRecHit = 
00262           (const SiTrackerGSMatchedRecHit2D*) (&(*theSeedItr));
00263         currentHit = TrackerRecHit(theSeedingRecHit,theTrackerGeometry);
00264         theHits.push_back(currentHit);
00265       }
00266       
00267       // Loop over clusters     
00268       unsigned csize = sclRefs.size();
00269       for  (unsigned int i=0;i<csize;++i) {
00270         
00271         // Find the pixel seeds (actually only the best one is returned)
00272 #ifdef FAMOS_DEBUG
00273         std::cout << "new cluster, calling addAseedFromThisCluster" << std::endl;
00274 #endif
00275         addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]);
00276         
00277       } 
00278       // End loop over clusters
00279     }
00280     // End loop over seeds
00281   }
00282   // end else
00283 
00284   // Back to the expected collection
00285  
00286   std::map<unsigned,std::vector<reco::ElectronPixelSeed> >::const_iterator px = myPixelSeeds.begin();
00287   std::map<unsigned,std::vector<reco::ElectronPixelSeed> >::const_iterator pxEnd = myPixelSeeds.end();
00288   for ( ; px!=pxEnd; ++px ) {
00289     unsigned nSeeds = (px->second).size();
00290     for ( unsigned ipx = 0; ipx<nSeeds; ++ipx ) { 
00291       out.push_back((px->second)[ipx]); 
00292       reco::ElectronPixelSeed is = px->second[ipx];
00293     }
00294   }
00295 
00296   LogDebug ("run") << ": For event "<<e.id();
00297   LogDebug ("run") <<"Nr of superclusters: "<<sclRefs.size()
00298                    <<", no. of ElectronPixelSeeds found  = " << out.size();
00299 #ifdef FAMOS_DEBUG
00300   std::cout << ": For event "<<e.id() << std::endl;
00301   std::cout <<"Nr of superclusters: "<<sclRefs.size()
00302             <<", no. of ElectronPixelSeeds found  = " << out.size() << std::endl;
00303 #endif
00304   
00305 }

void ElectronGSPixelSeedGenerator::setup ( bool   ) 

void ElectronGSPixelSeedGenerator::setupES ( const edm::EventSetup setup  ) 

Definition at line 109 of file ElectronGSPixelSeedGenerator.cc.

References alongMomentum, edm::EventSetup::get(), myGSPixelMatcher, GSPixelHitMatcher::setES(), theGeomSearchTracker, theMagField, theMagneticFieldMap, thePropagator, theSetup, theTrackerGeometry, and theTrackerInteractionGeometry.

Referenced by ElectronGSPixelSeedProducer::beginRun().

00109                                                                      {
00110 
00111   theSetup= &setup;
00112 
00113   edm::ESHandle<MagneticField> pMF;
00114   setup.get<IdealMagneticFieldRecord>().get(pMF);
00115   theMagField = &(*pMF);
00116 
00117   edm::ESHandle<TrackerGeometry>        geometry;
00118   setup.get<TrackerDigiGeometryRecord>().get(geometry);
00119   theTrackerGeometry = &(*geometry);
00120 
00121   edm::ESHandle<GeometricSearchTracker> recoGeom;
00122   setup.get<TrackerRecoGeometryRecord>().get( recoGeom );
00123   theGeomSearchTracker = &(*recoGeom);
00124 
00125   edm::ESHandle<TrackerInteractionGeometry> interGeom;
00126   setup.get<TrackerInteractionGeometryRecord>().get( interGeom );
00127   theTrackerInteractionGeometry = &(*interGeom);
00128 
00129   edm::ESHandle<MagneticFieldMap> fieldMap;
00130   setup.get<MagneticFieldMapRecord>().get(fieldMap);
00131   theMagneticFieldMap = &(*fieldMap);
00132 
00133   thePropagator = new PropagatorWithMaterial(alongMomentum,.000511,&(*theMagField)); 
00134 
00135   myGSPixelMatcher->setES(theMagneticFieldMap,
00136                           theTrackerGeometry,
00137                           theGeomSearchTracker,
00138                           theTrackerInteractionGeometry);
00139 
00140 }


Member Data Documentation

math::XYZPoint ElectronGSPixelSeedGenerator::BSPosition_ [private]

Definition at line 111 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster(), and run().

float ElectronGSPixelSeedGenerator::deltaPhi1High_ [private]

Definition at line 105 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

float ElectronGSPixelSeedGenerator::deltaPhi1Low_ [private]

Definition at line 105 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

float ElectronGSPixelSeedGenerator::deltaPhi2_ [private]

Definition at line 106 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

bool ElectronGSPixelSeedGenerator::dynamicphiroad_ [private]

Definition at line 99 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

bool ElectronGSPixelSeedGenerator::fromTrackerSeeds_ [private]

Definition at line 117 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster(), ElectronGSPixelSeedGenerator(), and run().

float ElectronGSPixelSeedGenerator::highPtThreshold_ [private]

Definition at line 102 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

float ElectronGSPixelSeedGenerator::lowPtThreshold_ [private]

Definition at line 101 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

GSPixelHitMatcher* ElectronGSPixelSeedGenerator::myGSPixelMatcher [private]

Definition at line 113 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster(), ElectronGSPixelSeedGenerator(), run(), setupES(), and ~ElectronGSPixelSeedGenerator().

float ElectronGSPixelSeedGenerator::phimax2_ [private]

Definition at line 104 of file ElectronGSPixelSeedGenerator.h.

float ElectronGSPixelSeedGenerator::phimin2_ [private]

Definition at line 104 of file ElectronGSPixelSeedGenerator.h.

double ElectronGSPixelSeedGenerator::pTMin2 [private]

Definition at line 110 of file ElectronGSPixelSeedGenerator.h.

Referenced by run().

PTrajectoryStateOnDet* ElectronGSPixelSeedGenerator::pts_ [private]

Definition at line 132 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster(), and prepareElTrackSeed().

PRecHitContainer ElectronGSPixelSeedGenerator::recHits_ [private]

Definition at line 131 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster(), and prepareElTrackSeed().

bool ElectronGSPixelSeedGenerator::searchInTIDTEC [private]

Definition at line 107 of file ElectronGSPixelSeedGenerator.h.

Referenced by ElectronGSPixelSeedGenerator(), and run().

float ElectronGSPixelSeedGenerator::sizeWindowENeg_ [private]

Definition at line 103 of file ElectronGSPixelSeedGenerator.h.

Referenced by addASeedToThisCluster().

const edm::InputTag ElectronGSPixelSeedGenerator::theBeamSpot [private]

Definition at line 129 of file ElectronGSPixelSeedGenerator.h.

Referenced by run().

const GeometricSearchTracker* ElectronGSPixelSeedGenerator::theGeomSearchTracker [private]

Definition at line 122 of file ElectronGSPixelSeedGenerator.h.

Referenced by setupES().

TrajectorySeedCollection* ElectronGSPixelSeedGenerator::theInitialSeedColl [private]

Definition at line 116 of file ElectronGSPixelSeedGenerator.h.

Referenced by run().

const MagneticField* ElectronGSPixelSeedGenerator::theMagField [private]

Definition at line 119 of file ElectronGSPixelSeedGenerator.h.

Referenced by setupES().

const MagneticFieldMap* ElectronGSPixelSeedGenerator::theMagneticFieldMap [private]

Definition at line 120 of file ElectronGSPixelSeedGenerator.h.

Referenced by setupES().

PropagatorWithMaterial* ElectronGSPixelSeedGenerator::thePropagator [private]

Definition at line 126 of file ElectronGSPixelSeedGenerator.h.

Referenced by prepareElTrackSeed(), setupES(), and ~ElectronGSPixelSeedGenerator().

const edm::EventSetup* ElectronGSPixelSeedGenerator::theSetup [private]

Definition at line 128 of file ElectronGSPixelSeedGenerator.h.

Referenced by prepareElTrackSeed(), and setupES().

const TrackerGeometry* ElectronGSPixelSeedGenerator::theTrackerGeometry [private]

Definition at line 121 of file ElectronGSPixelSeedGenerator.h.

Referenced by run(), and setupES().

const TrackerInteractionGeometry* ElectronGSPixelSeedGenerator::theTrackerInteractionGeometry [private]

Definition at line 123 of file ElectronGSPixelSeedGenerator.h.

Referenced by setupES().

KFUpdator* ElectronGSPixelSeedGenerator::theUpdator [private]

Definition at line 125 of file ElectronGSPixelSeedGenerator.h.

Referenced by ~ElectronGSPixelSeedGenerator().

TrajectoryStateTransform ElectronGSPixelSeedGenerator::transformer_ [private]

Definition at line 130 of file ElectronGSPixelSeedGenerator.h.

Referenced by prepareElTrackSeed().

double ElectronGSPixelSeedGenerator::zmax1_ [private]

Definition at line 109 of file ElectronGSPixelSeedGenerator.h.

double ElectronGSPixelSeedGenerator::zmin1_ [private]

Definition at line 109 of file ElectronGSPixelSeedGenerator.h.


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:20:13 2009 for CMSSW by  doxygen 1.5.4