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