00001
00002
00003
00004
00005
00011
00012
00013
00014 #include <vector>
00015 #include <utility>
00016
00017 #include "DataFormats/Math/interface/Point3D.h"
00018 #include "DataFormats/Common/interface/Handle.h"
00019
00020 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00021 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00022
00023 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00024
00025 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00026
00027 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00028
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00031
00032 #include "RecoTracker/TkSeedGenerator/interface/FastCircle.h"
00033 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripMatchedRecHit.h"
00034 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00035 #include "RecoTracker/TkSeedGenerator/interface/FastHelix.h"
00036 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00037
00038 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00039 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00040
00041 #include "MagneticField/Engine/interface/MagneticField.h"
00042 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00043
00044
00045
00046 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00047 #include "TrackingTools/DetLayers/interface/GeometricSearchDet.h"
00048 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00049 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00050 #include "TrackingTools/MeasurementDet/interface/MeasurementDet.h"
00051 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00052 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00053
00054 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00055
00056 #include "RecoEgamma/EgammaElectronAlgos/interface/SiStripElectronSeedGenerator.h"
00057
00058 SiStripElectronSeedGenerator::SiStripElectronSeedGenerator(const edm::ParameterSet &pset)
00059 : beamSpotTag_("offlineBeamSpot"),
00060 theUpdator(0),thePropagator(0),theMeasurementTracker(0),
00061 theSetup(0),pts_(0),theMatcher_(0),
00062 cacheIDMagField_(0),cacheIDCkfComp_(0),cacheIDTrkGeom_(0),
00063 tibOriginZCut_(pset.getParameter<double>("tibOriginZCut")),
00064 tidOriginZCut_(pset.getParameter<double>("tidOriginZCut")),
00065 tecOriginZCut_(pset.getParameter<double>("tecOriginZCut")),
00066 monoOriginZCut_(pset.getParameter<double>("monoOriginZCut")),
00067 tibDeltaPsiCut_(pset.getParameter<double>("tibDeltaPsiCut")),
00068 tidDeltaPsiCut_(pset.getParameter<double>("tidDeltaPsiCut")),
00069 tecDeltaPsiCut_(pset.getParameter<double>("tecDeltaPsiCut")),
00070 monoDeltaPsiCut_(pset.getParameter<double>("monoDeltaPsiCut")),
00071 tibPhiMissHit2Cut_(pset.getParameter<double>("tibPhiMissHit2Cut")),
00072 tidPhiMissHit2Cut_(pset.getParameter<double>("tidPhiMissHit2Cut")),
00073 tecPhiMissHit2Cut_(pset.getParameter<double>("tecPhiMissHit2Cut")),
00074 monoPhiMissHit2Cut_(pset.getParameter<double>("monoPhiMissHit2Cut")),
00075 tibZMissHit2Cut_(pset.getParameter<double>("tibZMissHit2Cut")),
00076 tidRMissHit2Cut_(pset.getParameter<double>("tidRMissHit2Cut")),
00077 tecRMissHit2Cut_(pset.getParameter<double>("tecRMissHit2Cut")),
00078 tidEtaUsage_(pset.getParameter<double>("tidEtaUsage")),
00079 tidMaxHits_(pset.getParameter<int>("tidMaxHits")),
00080 tecMaxHits_(pset.getParameter<int>("tecMaxHits")),
00081 monoMaxHits_(pset.getParameter<int>("monoMaxHits")),
00082 maxSeeds_(pset.getParameter<int>("maxSeeds"))
00083 {
00084
00085 if (pset.exists("measurementTrackerName"))
00086 { theMeasurementTrackerName = pset.getParameter<std::string>("measurementTrackerName") ; }
00087
00088
00089 if (pset.exists("beamSpot"))
00090 { beamSpotTag_ = pset.getParameter<edm::InputTag>("beamSpot") ; }
00091
00092 theUpdator = new KFUpdator();
00093 theEstimator = new Chi2MeasurementEstimator(30,3);
00094 }
00095
00096
00097 SiStripElectronSeedGenerator::~SiStripElectronSeedGenerator() {
00098 delete thePropagator;
00099 delete theUpdator;
00100 }
00101
00102
00103 void SiStripElectronSeedGenerator::setupES(const edm::EventSetup& setup) {
00104
00105 if (cacheIDMagField_!=setup.get<IdealMagneticFieldRecord>().cacheIdentifier()) {
00106 setup.get<IdealMagneticFieldRecord>().get(theMagField);
00107 cacheIDMagField_=setup.get<IdealMagneticFieldRecord>().cacheIdentifier();
00108 if (thePropagator) delete thePropagator;
00109 thePropagator = new PropagatorWithMaterial(alongMomentum,.000511,&(*theMagField));
00110 }
00111
00112 if (cacheIDCkfComp_!=setup.get<CkfComponentsRecord>().cacheIdentifier()) {
00113 setup.get<CkfComponentsRecord>().get(theMeasurementTrackerName,measurementTrackerHandle);
00114 cacheIDCkfComp_=setup.get<CkfComponentsRecord>().cacheIdentifier();
00115 theMeasurementTracker = measurementTrackerHandle.product();
00116 }
00117
00118 if (cacheIDTrkGeom_!=setup.get<TrackerDigiGeometryRecord>().cacheIdentifier()) {
00119 setup.get<TrackerDigiGeometryRecord>().get(trackerGeometryHandle);
00120 cacheIDTrkGeom_=setup.get<TrackerDigiGeometryRecord>().cacheIdentifier();
00121 }
00122
00123 }
00124
00125 void SiStripElectronSeedGenerator::run(edm::Event& e, const edm::EventSetup& setup,
00126 const edm::Handle<reco::SuperClusterCollection> &clusters,
00127 reco::ElectronSeedCollection & out) {
00128 theSetup= &setup;
00129 e.getByLabel(beamSpotTag_,theBeamSpot);
00130 theMeasurementTracker->update(e);
00131
00132 for (unsigned int i=0;i<clusters->size();++i) {
00133 edm::Ref<reco::SuperClusterCollection> theClusB(clusters,i);
00134
00135 LogDebug ("run") << "new cluster, calling findSeedsFromCluster";
00136 findSeedsFromCluster(theClusB,theBeamSpot,out);
00137 }
00138
00139 LogDebug ("run") << ": For event "<<e.id();
00140 LogDebug ("run") <<"Nr of superclusters: "<<clusters->size()
00141 <<", no. of ElectronSeeds found = " << out.size();
00142 }
00143
00144
00145
00146 void SiStripElectronSeedGenerator::findSeedsFromCluster
00147 ( edm::Ref<reco::SuperClusterCollection> seedCluster,
00148 edm::Handle<reco::BeamSpot> bs,
00149 reco::ElectronSeedCollection & result )
00150 {
00151
00152 layer1Hits_.clear() ;
00153 layer2Hits_.clear() ;
00154 backupLayer2Hits_.clear() ;
00155
00156 using namespace std;
00157
00158 double sCenergy = seedCluster->energy();
00159 math::XYZPoint sCposition = seedCluster->position();
00160 double scEta = seedCluster->eta();
00161
00162 double scz = sCposition.z();
00163 double scr = sqrt(pow(sCposition.x(),2)+pow(sCposition.y(),2));
00164
00165 double pT = sCenergy * seedCluster->position().rho()/sqrt(seedCluster->x()*seedCluster->x()+seedCluster->y()*seedCluster->y()+seedCluster->z()*seedCluster->z());
00166
00167 double magneticField = 3.8;
00168
00169
00170 double phiVsRSlope = -3.00e-3 * magneticField / pT / 2.;
00171
00172
00173
00174 GlobalPoint beamSpot(bs->x0(),bs->y0(),bs->z0());
00175 GlobalPoint superCluster(sCposition.x(),sCposition.y(),sCposition.z());
00176 double r0 = beamSpot.perp();
00177 double z0 = beamSpot.z();
00178
00179
00180
00181
00182
00183
00184 int chargeHypothesis;
00185 double chargeSelector = sCenergy - (int)sCenergy;
00186 if(chargeSelector >= 0.5) chargeHypothesis = -1;
00187 if(chargeSelector < 0.5) chargeHypothesis = 1;
00188
00189
00190 double rFake = 25.;
00191 double phiFake = phiDiff(superCluster.phi(),chargeHypothesis * phiVsRSlope * (scr - rFake));
00192 double zFake = (rFake*(scz-z0)-r0*scz+scr*z0)/(scr-r0);
00193 double xFake = rFake * cos(phiFake);
00194 double yFake = rFake * sin(phiFake);
00195 GlobalPoint fakePoint(xFake,yFake,zFake);
00196
00197
00198 FastHelix initialHelix(superCluster,fakePoint,beamSpot,*theSetup);
00199
00200
00201 FreeTrajectoryState initialFTS = initialHelix.stateAtVertex();
00202
00203
00204 TransverseImpactPointExtrapolator* tipe = new TransverseImpactPointExtrapolator(*thePropagator);
00205 TrajectoryStateOnSurface initialTSOS = tipe->extrapolate(initialFTS,beamSpot);
00206
00207
00208 const GeometricSearchTracker* gst = theMeasurementTracker->geometricSearchTracker();
00209
00210 std::vector<BarrelDetLayer*> tibLayers = gst->tibLayers();
00211 DetLayer* tib1 = tibLayers.at(0);
00212 DetLayer* tib2 = tibLayers.at(1);
00213
00214 std::vector<ForwardDetLayer*> tecLayers;
00215 std::vector<ForwardDetLayer*> tidLayers;
00216 if(scEta < 0){
00217 tecLayers = gst->negTecLayers();
00218 tidLayers = gst->negTidLayers();
00219 }
00220 if(scEta > 0){
00221 tecLayers = gst->posTecLayers();
00222 tidLayers = gst->posTidLayers();
00223 }
00224
00225 DetLayer* tid1 = tidLayers.at(0);
00226 DetLayer* tid2 = tidLayers.at(1);
00227 DetLayer* tid3 = tidLayers.at(2);
00228 DetLayer* tec1 = tecLayers.at(0);
00229 DetLayer* tec2 = tecLayers.at(1);
00230 DetLayer* tec3 = tecLayers.at(2);
00231
00232
00233 std::vector<bool> useDL = useDetLayer(scEta);
00234 bool useTID = false;
00235
00236
00237
00238 int tid1MHC = 0;
00239 int tid2MHC = 0;
00240 int tid3MHC = 0;
00241 int tid1BHC = 0;
00242 int tid2BHC = 0;
00243 int tid3BHC = 0;
00244 int tec1MHC = 0;
00245 int tec2MHC = 0;
00246 int tec3MHC = 0;
00247
00248
00249 int seedCounter = 0;
00250
00251 bool hasLay1Hit = false;
00252 bool hasLay2Hit = false;
00253 bool hasBackupHit = false;
00254
00255 LayerMeasurements layerMeasurements(theMeasurementTracker);
00256
00257 std::vector<TrajectoryMeasurement> tib1measurements;
00258 if(useDL.at(0)) tib1measurements = layerMeasurements.measurements(*tib1,initialTSOS,*thePropagator,*theEstimator);
00259 std::vector<TrajectoryMeasurement> tib2measurements;
00260 if(useDL.at(1)) tib2measurements = layerMeasurements.measurements(*tib2,initialTSOS,*thePropagator,*theEstimator);
00261
00262
00263
00264
00265
00266 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tib1measurements.begin(); tmIter != tib1measurements.end(); ++ tmIter){
00267 ConstRecHitPointer hit = tmIter->recHit();
00268 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00269 if(matchedHit){
00270 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00271 if(preselection(position, superCluster, phiVsRSlope, 1)){
00272 hasLay1Hit = true;
00273 layer1Hits_.push_back(matchedHit);
00274 }
00275 }
00276 }
00277
00278 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tib2measurements.begin(); tmIter != tib2measurements.end(); ++ tmIter){
00279 ConstRecHitPointer hit = tmIter->recHit();
00280 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00281 if(matchedHit){
00282 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00283 if(preselection(position, superCluster, phiVsRSlope, 1)){
00284 hasLay2Hit = true;
00285 layer2Hits_.push_back(matchedHit);
00286 }
00287 }
00288 }
00289
00290 if(!(hasLay1Hit && hasLay2Hit)) useTID = true;
00291 if(std::abs(scEta) > tidEtaUsage_) useTID = true;
00292 std::vector<TrajectoryMeasurement> tid1measurements;
00293 if(useDL.at(2) && useTID) tid1measurements = layerMeasurements.measurements(*tid1,initialTSOS,*thePropagator,*theEstimator);
00294 std::vector<TrajectoryMeasurement> tid2measurements;
00295 if(useDL.at(3) && useTID) tid2measurements = layerMeasurements.measurements(*tid2,initialTSOS,*thePropagator,*theEstimator);
00296 std::vector<TrajectoryMeasurement> tid3measurements;
00297 if(useDL.at(4) && useTID) tid3measurements = layerMeasurements.measurements(*tid3,initialTSOS,*thePropagator,*theEstimator);
00298
00299 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid1measurements.begin(); tmIter != tid1measurements.end(); ++ tmIter){
00300 if(tid1MHC < tidMaxHits_){
00301 ConstRecHitPointer hit = tmIter->recHit();
00302 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00303 if(matchedHit){
00304 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00305 if(preselection(position, superCluster, phiVsRSlope, 2)){
00306 tid1MHC++;
00307 hasLay1Hit = true;
00308 layer1Hits_.push_back(matchedHit);
00309 hasLay2Hit = true;
00310 layer2Hits_.push_back(matchedHit);
00311 }
00312 }else if(useDL.at(8) && tid1BHC < monoMaxHits_){
00313 const SiStripRecHit2D* backupHit = backupHitConverter(hit);
00314 if(backupHit){
00315 GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
00316 if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
00317 tid1BHC++;
00318 hasBackupHit = true;
00319 backupLayer2Hits_.push_back(backupHit);
00320 }
00321 }
00322 }
00323 }
00324 }
00325
00326 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid2measurements.begin(); tmIter != tid2measurements.end(); ++ tmIter){
00327 if(tid2MHC < tidMaxHits_){
00328 ConstRecHitPointer hit = tmIter->recHit();
00329 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00330 if(matchedHit){
00331 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00332 if(preselection(position, superCluster, phiVsRSlope, 2)){
00333 tid2MHC++;
00334 hasLay1Hit = true;
00335 layer1Hits_.push_back(matchedHit);
00336 hasLay2Hit = true;
00337 layer2Hits_.push_back(matchedHit);
00338 }
00339 }else if(useDL.at(8) && tid2BHC < monoMaxHits_){
00340 const SiStripRecHit2D* backupHit = backupHitConverter(hit);
00341 if(backupHit){
00342 GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
00343 if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
00344 tid2BHC++;
00345 hasBackupHit = true;
00346 backupLayer2Hits_.push_back(backupHit);
00347 }
00348 }
00349 }
00350 }
00351 }
00352
00353 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tid3measurements.begin(); tmIter != tid3measurements.end(); ++ tmIter){
00354 if(tid3MHC < tidMaxHits_){
00355 ConstRecHitPointer hit = tmIter->recHit();
00356 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00357 if(matchedHit){
00358 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00359 if(preselection(position, superCluster, phiVsRSlope, 2)){
00360 tid3MHC++;
00361 hasLay1Hit = true;
00362 layer1Hits_.push_back(matchedHit);
00363 hasLay2Hit = true;
00364 layer2Hits_.push_back(matchedHit);
00365 }
00366 }else if(useDL.at(8) && tid3BHC < monoMaxHits_){
00367 const SiStripRecHit2D* backupHit = backupHitConverter(hit);
00368 if(backupHit){
00369 GlobalPoint position = trackerGeometryHandle->idToDet(backupHit->geographicalId())->surface().toGlobal(backupHit->localPosition());
00370 if(preselection(position, superCluster, phiVsRSlope, 4) && position.perp() > 37.){
00371 tid3BHC++;
00372 hasBackupHit = true;
00373 backupLayer2Hits_.push_back(backupHit);
00374 }
00375 }
00376 }
00377 }
00378 }
00379
00380 std::vector<TrajectoryMeasurement> tec1measurements;
00381 if(useDL.at(5)) tec1measurements = layerMeasurements.measurements(*tec1,initialTSOS,*thePropagator,*theEstimator);
00382 std::vector<TrajectoryMeasurement> tec2measurements;
00383 if(useDL.at(6)) tec2measurements = layerMeasurements.measurements(*tec2,initialTSOS,*thePropagator,*theEstimator);
00384 std::vector<TrajectoryMeasurement> tec3measurements;
00385 if(useDL.at(7)) tec3measurements = layerMeasurements.measurements(*tec3,initialTSOS,*thePropagator,*theEstimator);
00386
00387 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec1measurements.begin(); tmIter != tec1measurements.end(); ++ tmIter){
00388 if(tec1MHC < tecMaxHits_){
00389 ConstRecHitPointer hit = tmIter->recHit();
00390 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00391 if(matchedHit){
00392 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00393 if(preselection(position, superCluster, phiVsRSlope, 3)){
00394 tec1MHC++;
00395 hasLay1Hit = true;
00396 layer1Hits_.push_back(matchedHit);
00397 hasLay2Hit = true;
00398 layer2Hits_.push_back(matchedHit);
00399 }
00400 }
00401 }
00402 }
00403
00404 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec2measurements.begin(); tmIter != tec2measurements.end(); ++ tmIter){
00405 if(tec2MHC < tecMaxHits_){
00406 ConstRecHitPointer hit = tmIter->recHit();
00407 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00408 if(matchedHit){
00409 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00410 if(preselection(position, superCluster, phiVsRSlope, 3)){
00411 tec2MHC++;
00412 hasLay1Hit = true;
00413 layer1Hits_.push_back(matchedHit);
00414 hasLay2Hit = true;
00415 layer2Hits_.push_back(matchedHit);
00416 }
00417 }
00418 }
00419 }
00420
00421 for(std::vector<TrajectoryMeasurement>::const_iterator tmIter = tec3measurements.begin(); tmIter != tec3measurements.end(); ++ tmIter){
00422 if(tec3MHC < tecMaxHits_){
00423 ConstRecHitPointer hit = tmIter->recHit();
00424 const SiStripMatchedRecHit2D* matchedHit = matchedHitConverter(hit);
00425 if(matchedHit){
00426 GlobalPoint position = trackerGeometryHandle->idToDet(matchedHit->geographicalId())->surface().toGlobal(matchedHit->localPosition());
00427 if(preselection(position, superCluster, phiVsRSlope, 3)){
00428 tec3MHC++;
00429 hasLay2Hit = true;
00430 layer2Hits_.push_back(matchedHit);
00431 }
00432 }
00433 }
00434 }
00435
00436
00437 if( hasLay1Hit && hasLay2Hit ){
00438
00439 for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1 = layer1Hits_.begin() ; hit1!= layer1Hits_.end(); ++hit1) {
00440 for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit2 = layer2Hits_.begin() ; hit2!= layer2Hits_.end(); ++hit2) {
00441
00442 if(seedCounter < maxSeeds_){
00443
00444 if(checkHitsAndTSOS(hit1,hit2,scr,scz,pT,scEta)) {
00445
00446 seedCounter++;
00447
00448 recHits_.clear();
00449
00450 SiStripMatchedRecHit2D *hit;
00451 hit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit1) ) );
00452 recHits_.push_back(hit);
00453 hit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit2) ) );
00454 recHits_.push_back(hit);
00455
00456 PropagationDirection dir = alongMomentum;
00457 reco::ElectronSeed seed(*pts_,recHits_,dir) ;
00458 reco::ElectronSeed::CaloClusterRef caloCluster(seedCluster) ;
00459 seed.setCaloCluster(caloCluster) ;
00460 result.push_back(seed);
00461
00462 delete pts_;
00463 pts_=0;
00464 }
00465
00466 }
00467
00468 }
00469
00470 }
00471
00472 }
00473
00474
00475
00476 if(hasLay1Hit && hasBackupHit && seedCounter == 0){
00477
00478 for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1 = layer1Hits_.begin() ; hit1!= layer1Hits_.end(); ++hit1) {
00479 for (std::vector<const SiStripRecHit2D*>::const_iterator hit2 = backupLayer2Hits_.begin() ; hit2!= backupLayer2Hits_.end(); ++hit2) {
00480
00481 if(seedCounter < maxSeeds_){
00482
00483 if(altCheckHitsAndTSOS(hit1,hit2,scr,scz,pT,scEta)) {
00484
00485 seedCounter++;
00486
00487 recHits_.clear();
00488
00489 SiStripMatchedRecHit2D *innerHit;
00490 innerHit=new SiStripMatchedRecHit2D(*(dynamic_cast <const SiStripMatchedRecHit2D *> (*hit1) ) );
00491 recHits_.push_back(innerHit);
00492 SiStripRecHit2D *outerHit;
00493 outerHit=new SiStripRecHit2D(*(dynamic_cast <const SiStripRecHit2D *> (*hit2) ) );
00494 recHits_.push_back(outerHit);
00495
00496 PropagationDirection dir = alongMomentum;
00497 reco::ElectronSeed seed(*pts_,recHits_,dir) ;
00498 reco::ElectronSeed::CaloClusterRef caloCluster(seedCluster) ;
00499 seed.setCaloCluster(caloCluster) ;
00500 result.push_back(seed);
00501
00502 delete pts_;
00503 pts_=0;
00504 }
00505
00506 }
00507
00508 }
00509
00510 }
00511
00512 }
00513
00514 }
00515
00516
00517
00518 bool SiStripElectronSeedGenerator::checkHitsAndTSOS(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1,
00519 std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit2,
00520 double rc,double zc,double pT,double scEta) {
00521
00522 bool seedCutsSatisfied = false;
00523
00524 using namespace std;
00525
00526 GlobalPoint hit1Pos = trackerGeometryHandle->idToDet((*hit1)->geographicalId())->surface().toGlobal((*hit1)->localPosition());
00527 double r1 = sqrt(hit1Pos.x()*hit1Pos.x() + hit1Pos.y()*hit1Pos.y());
00528 double phi1 = hit1Pos.phi();
00529 double z1=hit1Pos.z();
00530
00531 GlobalPoint hit2Pos = trackerGeometryHandle->idToDet((*hit2)->geographicalId())->surface().toGlobal((*hit2)->localPosition());
00532 double r2 = sqrt(hit2Pos.x()*hit2Pos.x() + hit2Pos.y()*hit2Pos.y());
00533 double phi2 = hit2Pos.phi();
00534 double z2 = hit2Pos.z();
00535
00536 if(r2 > r1 && (std::abs(z2) > std::abs(z1) || std::abs(scEta) < 0.25)) {
00537
00538
00539
00540 double curv = pT*100*.877;
00541
00542
00543 double a = (r2-r1)/(2*curv);
00544 double b = phiDiff(phi2,phi1);
00545
00546 double phiMissHit2=0;
00547 if(std::abs(b - a)<std::abs(b + a)) phiMissHit2 = b - a;
00548 if(std::abs(b - a)>std::abs(b + a)) phiMissHit2 = b + a;
00549
00550 double zMissHit2 = z2 - (r2*(zc-z1)-r1*zc+rc*z1)/(rc-r1);
00551
00552 double rPredHit2 = r1 + (rc-r1)/(zc-z1)*(z2-z1);
00553 double rMissHit2 = r2 - rPredHit2;
00554
00555 int subdetector = whichSubdetector(hit2);
00556
00557 bool zDiff = true;
00558 double zVar1 = std::abs(z1);
00559 double zVar2 = std::abs(z2 - z1);
00560 if(zVar1 > 75 && zVar1 < 95 && (zVar2 > 18 || zVar2 < 5)) zDiff = false;
00561 if(zVar1 > 100 && zVar1 < 110 && (zVar2 > 35 || zVar2 < 5)) zDiff = false;
00562 if(zVar1 > 125 && zVar1 < 150 && (zVar2 > 18 || zVar2 < 5)) zDiff = false;
00563
00564 if(subdetector == 1){
00565 int tibExtraCut = 0;
00566 if(r1 > 23 && r1 < 28 && r2 > 31 && r2 < 37) tibExtraCut = 1;
00567 if(std::abs(phiMissHit2) < tibPhiMissHit2Cut_ && std::abs(zMissHit2) < tibZMissHit2Cut_ && tibExtraCut == 1) seedCutsSatisfied = true;
00568 }else if(subdetector == 2){
00569 int tidExtraCut = 0;
00570 if(r1 > 23 && r1 < 34 && r2 > 26 && r2 < 42) tidExtraCut = 1;
00571 if(std::abs(phiMissHit2) < tidPhiMissHit2Cut_ && std::abs(rMissHit2) < tidRMissHit2Cut_ && tidExtraCut == 1 && zDiff) seedCutsSatisfied = true;
00572 }else if(subdetector == 3){
00573 int tecExtraCut = 0;
00574 if(r1 > 23 && r1 < 32 && r2 > 26 && r2 < 42) tecExtraCut = 1;
00575 if(std::abs(phiMissHit2) < tecPhiMissHit2Cut_ && std::abs(rMissHit2) < tecRMissHit2Cut_ && tecExtraCut == 1 && zDiff) seedCutsSatisfied = true;
00576 }
00577
00578 }
00579
00580 if(!seedCutsSatisfied) return false;
00581
00582
00583
00584 pts_=0;
00585
00586
00587
00588
00589
00590 RecHitPointer hit1Trans = TSiStripMatchedRecHit::build(trackerGeometryHandle->idToDet((*hit1)->geographicalId()), *hit1, theMatcher_);
00591 RecHitPointer hit2Trans = TSiStripMatchedRecHit::build(trackerGeometryHandle->idToDet((*hit2)->geographicalId()), *hit2, theMatcher_);
00592
00593 typedef TrajectoryStateOnSurface TSOS;
00594
00595 double vertexZ = z1 - (r1 * (zc - z1) ) / (rc - r1);
00596 GlobalPoint eleVertex(0.,0.,vertexZ);
00597
00598
00599 FastHelix helix(hit2Pos,hit1Pos,eleVertex,*theSetup);
00600 if (!helix.isValid()) return false;
00601
00602 FreeTrajectoryState fts = helix.stateAtVertex();
00603 TSOS propagatedState = thePropagator->propagate(fts,hit1Trans->det()->surface());
00604
00605 if (!propagatedState.isValid()) return false;
00606
00607 TSOS updatedState = theUpdator->update(propagatedState, *hit1Trans);
00608 TSOS propagatedState_out = thePropagator->propagate(fts,hit2Trans->det()->surface()) ;
00609
00610 if (!propagatedState_out.isValid()) return false;
00611
00612
00613
00614 TSOS updatedState_out = theUpdator->update(propagatedState_out, *hit2Trans);
00615
00616 pts_ = transformer_.persistentState(updatedState_out, hit2Trans->geographicalId().rawId());
00617
00618 return true;
00619 }
00620
00621 bool SiStripElectronSeedGenerator::altCheckHitsAndTSOS(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit1,
00622 std::vector<const SiStripRecHit2D*>::const_iterator hit2,
00623 double rc,double zc,double pT,double scEta) {
00624
00625 bool seedCutSatisfied = false;
00626
00627 using namespace std;
00628
00629 GlobalPoint hit1Pos = trackerGeometryHandle->idToDet((*hit1)->geographicalId())->surface().toGlobal((*hit1)->localPosition());
00630 double r1 = sqrt(hit1Pos.x()*hit1Pos.x() + hit1Pos.y()*hit1Pos.y());
00631 double phi1 = hit1Pos.phi();
00632 double z1=hit1Pos.z();
00633
00634 GlobalPoint hit2Pos = trackerGeometryHandle->idToDet((*hit2)->geographicalId())->surface().toGlobal((*hit2)->localPosition());
00635 double r2 = sqrt(hit2Pos.x()*hit2Pos.x() + hit2Pos.y()*hit2Pos.y());
00636 double phi2 = hit2Pos.phi();
00637 double z2 = hit2Pos.z();
00638
00639 if(r2 > r1 && std::abs(z2) > std::abs(z1)) {
00640
00641
00642
00643 double curv = pT*100*.877;
00644
00645
00646 double a = (r2-r1)/(2*curv);
00647 double b = phiDiff(phi2,phi1);
00648 double phiMissHit2 = 0;
00649 if(std::abs(b - a)<std::abs(b + a)) phiMissHit2 = b - a;
00650 if(std::abs(b - a)>std::abs(b + a)) phiMissHit2 = b + a;
00651
00652 if(std::abs(phiMissHit2) < monoPhiMissHit2Cut_) seedCutSatisfied = true;
00653
00654 }
00655
00656 if(!seedCutSatisfied) return false;
00657
00658
00659
00660 pts_=0;
00661
00662
00663
00664
00665
00666 RecHitPointer hit1Trans = TSiStripMatchedRecHit::build(trackerGeometryHandle->idToDet((*hit1)->geographicalId()), *hit1, theMatcher_);
00667 RecHitPointer hit2Trans = TSiStripMatchedRecHit::build(trackerGeometryHandle->idToDet((*hit2)->geographicalId()), *hit2, theMatcher_);
00668
00669 typedef TrajectoryStateOnSurface TSOS;
00670
00671 double vertexZ = z1 - (r1 * (zc - z1) ) / (rc - r1);
00672 GlobalPoint eleVertex(0.,0.,vertexZ);
00673
00674
00675 FastHelix helix(hit2Pos,hit1Pos,eleVertex,*theSetup);
00676 if (!helix.isValid()) return false;
00677
00678 FreeTrajectoryState fts = helix.stateAtVertex();
00679 TSOS propagatedState = thePropagator->propagate(fts,hit1Trans->det()->surface());
00680
00681 if (!propagatedState.isValid()) return false;
00682
00683 TSOS updatedState = theUpdator->update(propagatedState, *hit1Trans);
00684 TSOS propagatedState_out = thePropagator->propagate(fts,hit2Trans->det()->surface()) ;
00685
00686 if (!propagatedState_out.isValid()) return false;
00687
00688
00689
00690 TSOS updatedState_out = theUpdator->update(propagatedState_out, *hit2Trans);
00691
00692 pts_ = transformer_.persistentState(updatedState_out, hit2Trans->geographicalId().rawId());
00693
00694 return true;
00695 }
00696
00697
00698 bool SiStripElectronSeedGenerator::preselection(GlobalPoint position,GlobalPoint superCluster,double phiVsRSlope,int hitLayer){
00699 double r = position.perp();
00700 double phi = position.phi();
00701 double z = position.z();
00702 double scr = superCluster.perp();
00703 double scphi = superCluster.phi();
00704 double scz = superCluster.z();
00705 double psi = phiDiff(phi,scphi);
00706 double deltaPsi = psi - (scr-r)*phiVsRSlope;
00707 double antiDeltaPsi = psi - (r-scr)*phiVsRSlope;
00708 double dP;
00709 if (std::abs(deltaPsi)<std::abs(antiDeltaPsi)){
00710 dP = deltaPsi;
00711 }else{
00712 dP = antiDeltaPsi;
00713 }
00714 double originZ = (scr*z - r*scz)/(scr-r);
00715
00716 bool result = false;
00717
00718 if(hitLayer == 1){
00719 if(std::abs(originZ) < tibOriginZCut_ && std::abs(dP) < tibDeltaPsiCut_) result = true;
00720 }else if(hitLayer == 2){
00721 if(std::abs(originZ) < tidOriginZCut_ && std::abs(dP) < tidDeltaPsiCut_) result = true;
00722 }else if(hitLayer == 3){
00723 if(std::abs(originZ) < tecOriginZCut_ && std::abs(dP) < tecDeltaPsiCut_) result = true;
00724 }else if(hitLayer == 4){
00725 if(std::abs(originZ) < monoOriginZCut_ && std::abs(dP) < monoDeltaPsiCut_) result = true;
00726 }
00727
00728 return result;
00729 }
00730
00731
00732
00733 int SiStripElectronSeedGenerator::whichSubdetector(std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit){
00734 int result = 0;
00735 if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TIB){
00736 result = 1;
00737 }else if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TID){
00738 result = 2;
00739 }else if(((*hit)->geographicalId()).subdetId() == StripSubdetector::TEC){
00740 result = 3;
00741 }
00742 return result;
00743 }
00744
00745 const SiStripMatchedRecHit2D* SiStripElectronSeedGenerator::matchedHitConverter(ConstRecHitPointer crhp){
00746 const TrackingRecHit* trh = crhp->hit();
00747 const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(trh);
00748 return matchedHit;
00749 }
00750
00751 const SiStripRecHit2D* SiStripElectronSeedGenerator::backupHitConverter(ConstRecHitPointer crhp){
00752 const TrackingRecHit* trh = crhp->hit();
00753 const SiStripRecHit2D* backupHit = dynamic_cast<const SiStripRecHit2D*>(trh);
00754 return backupHit;
00755 }
00756
00757 std::vector<bool> SiStripElectronSeedGenerator::useDetLayer(double scEta){
00758 std::vector<bool> useDetLayer;
00759 double variable = std::abs(scEta);
00760 if(variable > 0 && variable < 1.8){
00761 useDetLayer.push_back(true);
00762 }else{
00763 useDetLayer.push_back(false);
00764 }
00765 if(variable > 0 && variable < 1.5){
00766 useDetLayer.push_back(true);
00767 }else{
00768 useDetLayer.push_back(false);
00769 }
00770 if(variable > 1 && variable < 2.1){
00771 useDetLayer.push_back(true);
00772 }else{
00773 useDetLayer.push_back(false);
00774 }
00775 if(variable > 1 && variable < 2.2){
00776 useDetLayer.push_back(true);
00777 }else{
00778 useDetLayer.push_back(false);
00779 }
00780 if(variable > 1 && variable < 2.3){
00781 useDetLayer.push_back(true);
00782 }else{
00783 useDetLayer.push_back(false);
00784 }
00785 if(variable > 1.8 && variable < 2.5){
00786 useDetLayer.push_back(true);
00787 }else{
00788 useDetLayer.push_back(false);
00789 }
00790 if(variable > 1.8 && variable < 2.5){
00791 useDetLayer.push_back(true);
00792 }else{
00793 useDetLayer.push_back(false);
00794 }
00795 if(variable > 1.8 && variable < 2.5){
00796 useDetLayer.push_back(true);
00797 }else{
00798 useDetLayer.push_back(false);
00799 }
00800 if(variable > 1.2 && variable < 1.6){
00801 useDetLayer.push_back(true);
00802 }else{
00803 useDetLayer.push_back(false);
00804 }
00805 return useDetLayer;
00806 }
00807
00808
00809