CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoParticleFlow/PFSimProducer/plugins/ConvBremSeedProducer.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFSimProducer/plugins/ConvBremSeedProducer.h"
00002 #include "FWCore/Framework/interface/ESHandle.h"
00003 
00005 #include "FastSimulation/TrackerSetup/interface/TrackerInteractionGeometryRecord.h"
00006 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00007 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00008 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00009 #include "FastSimulation/ParticlePropagator/interface/MagneticFieldMapRecord.h"
00010 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00011 
00013 #include "FastSimulation/TrackerSetup/interface/TrackerInteractionGeometry.h"
00014 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00015 #include "MagneticField/Engine/interface/MagneticField.h"
00016 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
00017 
00019 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00020 #include "DataFormats/ParticleFlowReco/interface/ConvBremSeed.h"
00021 #include "DataFormats/ParticleFlowReco/interface/ConvBremSeedFwd.h"
00022 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrackFwd.h"
00023 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrack.h"
00024 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00025 #include "DataFormats/ParticleFlowReco/interface/PFBrem.h"
00026 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00027 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00028 
00030 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00031 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00032 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00033 #include "DataFormats/GeometrySurface/interface/Surface.h"
00034 #include "DataFormats/GeometrySurface/interface/TangentPlane.h"
00035 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
00036 #include "FastSimulation/ParticlePropagator/src/ParticlePropagator.cc"
00037 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00038 #include "FastSimulation/TrajectoryManager/interface/InsideBoundsMeasurementEstimator.h"
00039 #include "FastSimulation/TrajectoryManager/interface/LocalMagneticField.h"
00040 #include "RecoTracker/TkSeedGenerator/interface/FastHelix.h"
00041 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00042 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00043 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00044 
00045 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00046 
00047 using namespace edm;
00048 using namespace std;
00049 using namespace reco;
00050 
00051 ConvBremSeedProducer::ConvBremSeedProducer(const ParameterSet& iConfig):
00052   conf_(iConfig),
00053   fieldMap_(0),
00054   layerMap_(56, static_cast<const DetLayer*>(0)),
00055   negLayerOffset_(27)
00056 {
00057   produces<ConvBremSeedCollection>();
00058 }
00059 
00060 
00061 ConvBremSeedProducer::~ConvBremSeedProducer()
00062 {
00063  
00064 
00065 }
00066 
00067 
00068 void
00069 ConvBremSeedProducer::produce(Event& iEvent, const EventSetup& iSetup)
00070 {
00071  
00072   LogDebug("ConvBremSeedProducerProducer")<<"START event: "<<iEvent.id().event()
00073                                           <<" in run "<<iEvent.id().run();
00074 
00075   float pfmass=  0.0005;
00076 
00078 
00080   Handle<PFClusterCollection> PfC;
00081   iEvent.getByLabel(conf_.getParameter<InputTag>("PFClusters") ,PfC);
00082   const PFClusterCollection& PPP= *(PfC.product());
00083  
00085   Handle<SiPixelRecHitCollection> pixelHits;
00086   iEvent.getByLabel(conf_.getParameter<InputTag>("pixelRecHits") , pixelHits);
00087 
00089   Handle<SiStripRecHit2DCollection> rphirecHits;
00090   iEvent.getByLabel(conf_.getParameter<InputTag>("rphirecHits"),rphirecHits);
00091   Handle<SiStripRecHit2DCollection> stereorecHits;
00092   iEvent.getByLabel(conf_.getParameter<InputTag>("stereorecHits"), stereorecHits);
00093   Handle<SiStripMatchedRecHit2DCollection> matchedrecHits;
00094   iEvent.getByLabel(conf_.getParameter<InputTag>("matchedrecHits"), matchedrecHits);
00095 
00096   //GSFPFRECTRACKS
00097   Handle<GsfPFRecTrackCollection> thePfRecTrackCollection;
00098   iEvent.getByLabel(conf_.getParameter<InputTag>("PFRecTrackLabel"),
00099                     thePfRecTrackCollection);
00100   const GsfPFRecTrackCollection PfRTkColl = *(thePfRecTrackCollection.product()); 
00101 
00102 
00104   std::auto_ptr<ConvBremSeedCollection> output(new ConvBremSeedCollection);
00105 
00106 
00108   vector<pair< TrajectorySeed , pair<GlobalVector,float> > > unclean;
00109   //TRIPLET OF MODULES TO BE USED FOR SEEDING
00110   vector< vector< long int > > tripl;
00111   //LAYER MAP
00112   initializeLayerMap(); 
00113 
00114 
00115 
00117 
00118   for(unsigned int ipft=0;ipft<PfRTkColl.size();ipft++){
00119     GsfPFRecTrackRef pft(thePfRecTrackCollection,ipft);
00120     LogDebug("ConvBremSeedProducerProducer")<<"NEW GsfPFRecTRACK ";
00121     float eta_br=0;
00122     unclean.clear();
00123     tripl.clear();
00124     vector<int> gc;
00125     TrackingRecHitRefVector gsfRecHits=pft->gsfTrackRef()->extra()->recHits();
00126     float pfoutenergy=sqrt((pfmass*pfmass)+pft->gsfTrackRef()->outerMomentum().Mag2());
00127     XYZTLorentzVector mom =XYZTLorentzVector(pft->gsfTrackRef()->outerMomentum().x(),
00128                                              pft->gsfTrackRef()->outerMomentum().y(),
00129                                              pft->gsfTrackRef()->outerMomentum().z(),
00130                                              pfoutenergy);
00131     XYZTLorentzVector pos =   XYZTLorentzVector(pft->gsfTrackRef()->outerPosition().x(),
00132                                                 pft->gsfTrackRef()->outerPosition().y(),
00133                                                 pft->gsfTrackRef()->outerPosition().z(),
00134                                                 0.);
00135     BaseParticlePropagator theOutParticle=BaseParticlePropagator( RawParticle(mom,pos),
00136                                                                   0,0,B_.z());
00137     theOutParticle.setCharge(pft->gsfTrackRef()->charge());
00138 
00140     gc.push_back(GoodCluster(theOutParticle,PPP,0.5));
00141 
00142 
00143     vector<PFBrem> brem =(*pft).PFRecBrem();
00144     vector<PFBrem>::iterator ib=brem.begin();
00145     vector<PFBrem>::iterator ib_end=brem.end();
00146     LogDebug("ConvBremSeedProducerProducer")<<"NUMBER OF BREMS "<<brem.size();
00147 
00149     for (;ib!=ib_end;++ib){
00150     
00151       XYZTLorentzVector mom=pft->trajectoryPoint(ib->indTrajPoint()).momentum();
00152       XYZTLorentzVector pos=
00153         XYZTLorentzVector(
00154                           pft->trajectoryPoint(ib->indTrajPoint()).position().x(),
00155                           pft->trajectoryPoint(ib->indTrajPoint()).position().y(),
00156                           pft->trajectoryPoint(ib->indTrajPoint()).position().z(),
00157                           0);
00158       
00160       if (pos.Rho()>80) continue;
00161       if ((pos.Rho()>5)&&(fabs(ib->SigmaDeltaP()/ib->DeltaP())>3)) continue;  
00162       if (fabs(ib->DeltaP())<3) continue;
00163       eta_br=mom.eta();
00164       vector< vector< long int > >Idd;
00165       
00166 
00167 
00168       BaseParticlePropagator p( RawParticle(mom,pos),
00169                                 0,0,B_.z());
00170       p.setCharge(0);
00171       gc.push_back(GoodCluster(p,PPP,0.2));
00172 
00173       ParticlePropagator PP(p,fieldMap_);
00174 
00176       list<TrackerLayer>::const_iterator cyliter= geometry_->cylinderBegin();
00177       for ( ; cyliter != geometry_->cylinderEnd() ; ++cyliter ) {
00178 
00180        if (!(cyliter->sensitive())) continue;
00181        PP.setPropagationConditions(*cyliter);
00182        PP.propagate();    
00183        if (PP.getSuccess()==0) continue;
00184 
00186        LocalMagneticField mf(PP.getMagneticField());
00187        AnalyticalPropagator alongProp(&mf, anyDirection);
00188        InsideBoundsMeasurementEstimator est;
00189        const DetLayer* tkLayer = detLayer(*cyliter,PP.Z());
00190        if (&(*tkLayer)==0) continue;
00191        TrajectoryStateOnSurface trajState = makeTrajectoryState( tkLayer, PP, &mf);
00192             
00193        std::vector<DetWithState> compat 
00194          = tkLayer->compatibleDets( trajState, alongProp, est);
00195        vector <long int> temp;
00196        if (compat.size()==0) continue;
00197 
00198        for (std::vector<DetWithState>::const_iterator i=compat.begin(); i!=compat.end(); i++) {
00199               
00200          long int detid=i->first->geographicalId().rawId();
00201          
00202          if ((tkLayer->subDetector()!=GeomDetEnumerators::PixelBarrel)&&
00203              (tkLayer->subDetector()!=GeomDetEnumerators::PixelEndcap)){
00204            
00205             
00206            StDetMatch DetMatch = (rphirecHits.product())->find((detid));
00207            MatDetMatch MDetMatch =(matchedrecHits.product())->find((detid));
00208            
00209                 
00210            long int DetID=(DetMatch != rphirecHits->end())? detid:0;
00211 
00212            if ((MDetMatch != matchedrecHits->end()) && !MDetMatch->empty()) {
00213              long int pii = MDetMatch->begin()->monoId();
00214              StDetMatch CDetMatch = (rphirecHits.product())->find((pii));
00215              DetID=(CDetMatch != rphirecHits->end())? pii:0;
00216              
00217            }
00218 
00219            temp.push_back(DetID);
00220            
00221          }
00222          else{
00223            PiDetMatch DetMatch = (pixelHits.product())->find((detid));
00224            long int DetID=(DetMatch != pixelHits->end())? detid:0;
00225            temp.push_back(DetID);
00226         
00227         
00228          }
00229        }
00230 
00231        Idd.push_back(temp);
00232 
00233       }//END TRACKER LAYER LOOP
00234       if(Idd.size()<2)continue;
00235 
00237       for (unsigned int i=0;i<Idd.size()-2;i++){
00238         for (unsigned int i1=0;i1<Idd[i].size();i1++){
00239           for (unsigned int i2=0;i2<Idd[i+1].size();i2++){
00240             for (unsigned int i3=0;i3<Idd[i+2].size();i3++){
00241               if ((Idd[i][i1]!=0) &&(Idd[i+1][i2]!=0) &&(Idd[i+2][i3]!=0) ){
00242                 vector<long int >tmp;
00243                 tmp.push_back(Idd[i][i1]);  tmp.push_back(Idd[i+1][i2]); tmp.push_back(Idd[i+2][i3]);
00244                 
00245                 bool newTrip=true;
00246                 for (unsigned int iv=0;iv<tripl.size();iv++){
00247                   if((tripl[iv][0]==tmp[0])&&(tripl[iv][1]==tmp[1])&&(tripl[iv][2]==tmp[2])) newTrip=false;
00248 
00249                 }
00250                 if (newTrip){
00251 
00252                   tripl.push_back(tmp);
00253                 }
00254               }
00255             }     
00256           }       
00257         }
00258       }
00259     }//END BREM LOOP
00260 
00261     float sineta_brem =sinh(eta_br);
00262     
00263 
00264     //OUTPUT COLLECTION
00265 
00266     TransientTrackingRecHit::ConstRecHitContainer glob_hits;
00267     OwnVector<TrackingRecHit> loc_hits;
00268     for (unsigned int i=0;i<tripl.size();i++){
00269       StDetMatch DetMatch1 = (rphirecHits.product())->find(tripl[i][0]);
00270       StDetMatch DetMatch2 = (rphirecHits.product())->find(tripl[i][1]);
00271       StDetMatch DetMatch3 = (rphirecHits.product())->find(tripl[i][2]);
00272       if ((DetMatch1 == rphirecHits->end()) ||
00273           (DetMatch2 == rphirecHits->end()) ||
00274           (DetMatch3 == rphirecHits->end()) )  continue;
00275       StDetSet DetSet1 = *DetMatch1;
00276       StDetSet DetSet2 = *DetMatch2;
00277       StDetSet DetSet3 = *DetMatch3;
00278 
00279       for (StDetSet::const_iterator it1=DetSet1.begin();it1!=DetSet1.end();++it1){
00280         GlobalPoint gp1=tracker_->idToDet(tripl[i][0])->surface().
00281           toGlobal(it1->localPosition());
00282 
00283         bool tak1=isGsfTrack(gsfRecHits,&(*it1));
00284         
00285         for (StDetSet::const_iterator it2=DetSet2.begin();it2!=DetSet2.end();++it2){
00286           GlobalPoint gp2=tracker_->idToDet(tripl[i][1])->surface().
00287             toGlobal(it2->localPosition());
00288           bool tak2=isGsfTrack(gsfRecHits,&(*it2));
00289 
00290           for (StDetSet::const_iterator it3=DetSet3.begin();it3!=DetSet3.end();++it3){
00291             //  ips++;
00292             GlobalPoint gp3=tracker_->idToDet(tripl[i][2])->surface().
00293               toGlobal(it3->localPosition());    
00294             bool tak3=isGsfTrack(gsfRecHits,&(*it3));  
00295             
00296 
00297             FastHelix helix(gp3, gp2, gp1,iSetup);
00298             GlobalVector gv=helix.stateAtVertex().parameters().momentum();
00299             GlobalVector gv_corr(gv.x(),gv.y(),gv.perp()*sineta_brem);
00300             float ene= sqrt(gv_corr.mag2()+(pfmass*pfmass));
00301 
00302             GlobalPoint gp=helix.stateAtVertex().parameters().position();
00303             float ch=helix.stateAtVertex().parameters().charge();
00304 
00305 
00306 
00307 
00308             XYZTLorentzVector mom = XYZTLorentzVector(gv.x(),gv.y(),gv_corr.z(),ene);
00309             XYZTLorentzVector pos = XYZTLorentzVector(gp.x(),gp.y(),gp.z(),0.);
00310             BaseParticlePropagator theOutParticle(RawParticle(mom,pos),0,0,B_.z());
00311             theOutParticle.setCharge(ch);
00312             int bgc=GoodCluster(theOutParticle,PPP,0.3,true);
00313 
00314             if (gv.perp()<0.5) continue;
00315 
00316             if (tak1+tak2+tak3>2) continue;
00317 
00318             if (bgc==-1) continue;
00319             bool clTak=false;
00320             for (unsigned int igcc=0; igcc<gc.size(); igcc++){
00321               if (clTak) continue;
00322               if (bgc==gc[igcc]) clTak=true;
00323             }
00324             if (clTak) continue;
00325 
00326 
00327 
00328 
00329             GlobalTrajectoryParameters Gtp(gp1,
00330                                            gv,int(ch), 
00331                                            &(*magfield_));
00332             glob_hits.clear(); loc_hits.clear();
00333             glob_hits.push_back(hitBuilder_->build(it1->clone()));
00334             glob_hits.push_back(hitBuilder_->build(it2->clone()));
00335             glob_hits.push_back(hitBuilder_->build(it3->clone()));
00336 
00338 
00339             FreeTrajectoryState CSeed(Gtp,
00340                                       CurvilinearTrajectoryError(AlgebraicSymMatrix55(AlgebraicMatrixID())));
00341             TrajectoryStateOnSurface updatedState;        
00342   
00343             for (int ih=0;ih<3;ih++){
00344 
00345               TrajectoryStateOnSurface state = (ih==0)?
00346                 propagator_->propagate(CSeed,
00347                                        tracker_->idToDet(tripl[i][ih])->surface()):
00348                 propagator_->propagate(updatedState,
00349                                        tracker_->idToDet(tripl[i][ih])->surface());
00350               
00351               if (!state.isValid()){ 
00352                 ih=3;
00353                 continue;}
00354 
00355               updatedState =  kfUpdator_->update(state, *glob_hits[ih]);
00356               loc_hits.push_back(glob_hits[ih]->hit()->clone());
00357               if (ih==2){
00358                 PTrajectoryStateOnDet const & PTraj= 
00359                   trajectoryStateTransform::persistentState(updatedState,tripl[i][2]);
00360                 //              output->push_back(Trajectoryseed(PTraj,loc_hits,alongMomentum));
00361                 unclean.push_back(make_pair(TrajectorySeed(PTraj,loc_hits,alongMomentum), 
00362                                             make_pair(gv_corr,ch)));
00363               }
00364               //    }
00365 
00366             }
00367           }    
00368         }    
00369       }
00370     }
00371     vector<bool> inPhot = sharedHits(unclean);
00372     for (unsigned int iu=0; iu<unclean.size();iu++){
00373 
00374       if (inPhot[iu])
00375         output->push_back(ConvBremSeed(unclean[iu].first,pft));
00376 
00377     }
00378 
00379   } //END GSF TRACK COLLECTION LOOP 
00380   LogDebug("ConvBremSeedProducerProducer")<<"END";
00381   iEvent.put(output);
00382     
00383 }
00384 
00385 
00386 void 
00387 ConvBremSeedProducer::beginRun(edm::Run& run,
00388                                const EventSetup& iSetup)
00389 {
00390   ESHandle<GeometricSearchTracker> track;
00391   iSetup.get<TrackerRecoGeometryRecord>().get( track ); 
00392   geomSearchTracker_ = track.product();
00393   
00394   ESHandle<TrackerInteractionGeometry>  theTrackerInteractionGeometry;
00395   iSetup.get<TrackerInteractionGeometryRecord>().get(theTrackerInteractionGeometry );
00396   geometry_=theTrackerInteractionGeometry.product();
00397 
00398   ESHandle<TrackerGeometry> tracker;
00399   iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
00400   tracker_=tracker.product();
00401 
00402   ESHandle<MagneticField> magfield;
00403   iSetup.get<IdealMagneticFieldRecord>().get(magfield);
00404   magfield_=magfield.product();
00405   B_=magfield_->inTesla(GlobalPoint(0,0,0));
00406 
00407   ESHandle<MagneticFieldMap> fieldMap;
00408   iSetup.get<MagneticFieldMapRecord>().get(fieldMap);
00409   fieldMap_ =fieldMap.product();
00410 
00411 
00412 
00413   ESHandle<TransientTrackingRecHitBuilder> hitBuilder;
00414   iSetup.get<TransientRecHitRecord>().get(conf_.getParameter<string>("TTRHBuilder"),hitBuilder);
00415   hitBuilder_= hitBuilder.product();
00416 
00417   propagator_  =    new PropagatorWithMaterial(alongMomentum,0.0005,&(*magfield) );
00418   kfUpdator_   =    new KFUpdator();
00419 }
00420 
00421 void 
00422 ConvBremSeedProducer::endRun() {
00423   delete propagator_;
00424   delete kfUpdator_;
00425 }
00426 
00427 void 
00428 ConvBremSeedProducer::initializeLayerMap()
00429 {
00430 
00431 
00432   // These are the BoundSurface&, the BoundDisk* and the BoundCylinder* for that layer
00433   //   const BoundSurface& theSurface = layer.surface();
00434   //   BoundDisk* theDisk = layer.disk();  // non zero for endcaps
00435   //   BoundCylinder* theCylinder = layer.cylinder(); // non zero for barrel
00436   //   int theLayer = layer.layerNumber(); // 1->3 PixB, 4->5 PixD, 
00437   //                                       // 6->9 TIB, 10->12 TID, 
00438   //                                       // 13->18 TOB, 19->27 TEC
00439   
00442     
00443     std::vector< BarrelDetLayer*>   barrelLayers = 
00444       geomSearchTracker_->barrelLayers();
00445     LogDebug("FastTracker") << "Barrel DetLayer dump: ";
00446     for (std::vector< BarrelDetLayer*>::const_iterator bl=barrelLayers.begin();
00447          bl != barrelLayers.end(); ++bl) {
00448       LogDebug("FastTracker")<< "radius " << (**bl).specificSurface().radius(); 
00449     }
00450 
00451   std::vector< ForwardDetLayer*>  posForwardLayers = 
00452     geomSearchTracker_->posForwardLayers();
00453   LogDebug("FastTracker") << "Positive Forward DetLayer dump: ";
00454   for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00455        fl != posForwardLayers.end(); ++fl) {
00456     LogDebug("FastTracker") << "Z pos "
00457                             << (**fl).surface().position().z()
00458                             << " radii " 
00459                             << (**fl).specificSurface().innerRadius() 
00460                             << ", " 
00461                             << (**fl).specificSurface().outerRadius(); 
00462   }
00463 
00464   const float rTolerance = 1.5;
00465   const float zTolerance = 3.;
00466 
00467   LogDebug("FastTracker")<< "Dump of TrackerInteractionGeometry cylinders:";
00468   for( std::list<TrackerLayer>::const_iterator i=geometry_->cylinderBegin();
00469        i!=geometry_->cylinderEnd(); ++i) {
00470     const BoundCylinder* cyl = i->cylinder();
00471     const BoundDisk* disk = i->disk();
00472 
00473     LogDebug("FastTracker") << "Famos Layer no " << i->layerNumber()
00474                             << " is sensitive? " << i->sensitive()
00475                             << " pos " << i->surface().position();
00476     if (!i->sensitive()) continue;
00477 
00478     if (cyl != 0) {
00479 
00480       LogDebug("FastTracker") << " cylinder radius " << cyl->radius();
00481       bool found = false;
00482 
00483       for (std::vector< BarrelDetLayer*>::const_iterator 
00484              bl=barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
00485 
00486         if (fabs( cyl->radius() - (**bl).specificSurface().radius()) < rTolerance) {
00487 
00488           layerMap_[i->layerNumber()] = *bl;
00489           found = true;
00490           LogDebug("FastTracker")<< "Corresponding DetLayer found with radius "
00491                                  << (**bl).specificSurface().radius();
00492                   
00493           break;
00494         }
00495       }
00496       if (!found) {
00497         LogError("FastTracker") << "FAILED to find a corresponding DetLayer!";
00498       }
00499     }
00500     else {
00501       LogDebug("FastTracker") << " disk radii " << disk->innerRadius() 
00502                  << ", " << disk->outerRadius();
00503 
00504       bool found = false;
00505 
00506       for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00507            fl != posForwardLayers.end(); ++fl) {
00508         if (fabs( disk->position().z() - (**fl).surface().position().z()) < zTolerance) {
00509           layerMap_[i->layerNumber()] = *fl;
00510           found = true;
00511           LogDebug("FastTracker") << "Corresponding DetLayer found with Z pos "
00512                                   << (**fl).surface().position().z()
00513                                   << " and radii " 
00514                                   << (**fl).specificSurface().innerRadius() 
00515                                   << ", " 
00516                                   << (**fl).specificSurface().outerRadius(); 
00517           break;
00518         }
00519       }
00520       if (!found) {
00521         LogError("FastTracker") << "FAILED to find a corresponding DetLayer!";
00522       }
00523     }
00524   }
00525 
00526 }
00527 const DetLayer*  ConvBremSeedProducer::detLayer( const TrackerLayer& layer, float zpos) const
00528 {
00529   if (zpos > 0 || !layer.forward() ) return layerMap_[layer.layerNumber()];
00530   else return layerMap_[layer.layerNumber()+negLayerOffset_];
00531 }
00532 
00533 TrajectoryStateOnSurface 
00534 ConvBremSeedProducer::makeTrajectoryState( const DetLayer* layer, 
00535                                const ParticlePropagator& pp,
00536                                const MagneticField* field) const
00537 {
00538 
00539   GlobalPoint  pos( pp.X(), pp.Y(), pp.Z());
00540   GlobalVector mom( pp.Px(), pp.Py(), pp.Pz());
00541 
00542   ReferenceCountingPointer<TangentPlane> plane = layer->surface().tangentPlane(pos);
00543 
00544   return TrajectoryStateOnSurface
00545     (GlobalTrajectoryParameters( pos, mom, TrackCharge( pp.charge()), field), *plane);
00546 }
00547 bool ConvBremSeedProducer::isGsfTrack(TrackingRecHitRefVector tkv, const TrackingRecHit *h ){
00548   trackingRecHit_iterator ib=tkv.begin();
00549   trackingRecHit_iterator ie=tkv.end();
00550   bool istaken=false;
00551   //  for (;ib!=ie-2;++ib){
00552     for (;ib!=ie;++ib){
00553     if (istaken) continue;
00554     if (!((*ib)->isValid())) continue;
00555  
00556     istaken = (*ib)->sharesInput(h,TrackingRecHit::all);
00557   }
00558   return istaken;
00559 }
00560 vector<bool> ConvBremSeedProducer::sharedHits( vector<pair< TrajectorySeed, 
00561                                    pair<GlobalVector,float> > > unclean){
00562 
00563   vector<bool> goodseed;
00564   goodseed.clear();
00565   if (unclean.size()<2){
00566     for (unsigned int i=0;i<unclean.size();i++)
00567       goodseed.push_back(true);
00568   }else{
00569  
00570     for (unsigned int i=0;i<unclean.size();i++)
00571       goodseed.push_back(true);
00572     
00573     for (unsigned int iu=0; iu<unclean.size()-1;iu++){   
00574       if (!goodseed[iu]) continue;
00575       for (unsigned int iu2=iu+1; iu2<unclean.size();iu2++){
00576         if (!goodseed[iu]) continue;
00577         if (!goodseed[iu2]) continue;
00578       //    if (unclean[iu].second.second *unclean[iu2].second.second >0)continue;
00579         
00580         TrajectorySeed::const_iterator sh = unclean[iu].first.recHits().first;
00581         TrajectorySeed::const_iterator sh_end = unclean[iu].first.recHits().second; 
00582         
00583         unsigned int shar =0;
00584         for (;sh!=sh_end;++sh){ 
00585 
00586           TrajectorySeed::const_iterator sh2 = unclean[iu2].first.recHits().first;
00587           TrajectorySeed::const_iterator sh2_end = unclean[iu2].first.recHits().second; 
00588         for (;sh2!=sh2_end;++sh2){
00589 
00590           if ((*sh).sharesInput(&(*sh2),TrackingRecHit::all))
00591 
00592             shar++;
00593 
00594         }
00595         }
00596         if (shar>=2){
00597         if (unclean[iu].second.first.perp()<unclean[iu2].second.first.perp()) goodseed[iu]=false;
00598         else goodseed[iu2]=false;
00599         }
00600       
00601       } 
00602  
00603     }
00604   }
00605   return goodseed;
00606 }
00607 
00608 
00609 
00610 int ConvBremSeedProducer::GoodCluster(BaseParticlePropagator bpg, const PFClusterCollection& pfc, float minep, bool sec){
00611 
00612   bpg.propagateToEcalEntrance(false);
00613   float dr=1000;
00614   float de=1000;
00615   float df=1000;
00616   int ibest=-1;
00617 
00618   if(bpg.getSuccess()!=0){
00619 
00620     for (unsigned int i =0; i<pfc.size();i++ ){
00621       float tmp_ep=pfc[i].energy()/bpg.momentum().e();
00622       float tmp_phi=fabs(pfc[i].position().phi()-bpg.vertex().phi());
00623       if (tmp_phi>TMath::TwoPi()) tmp_phi-= TMath::TwoPi(); 
00624       float tmp_eta=fabs(pfc[i].position().eta()-bpg.vertex().eta());
00625       float tmp_dr=sqrt(pow(tmp_phi,2)+pow(tmp_eta,2));
00626       bool isBet=(tmp_dr<dr);
00627       if (sec) isBet=(tmp_phi<df);
00628       if ((isBet)&&(tmp_ep>minep)&&(tmp_ep<10)){
00629         dr=tmp_dr;
00630         de=tmp_eta;
00631         df=tmp_phi; 
00632         ibest=i;
00633       }
00634     }
00635     bool isBad=(dr>0.1);
00636     if (sec) isBad= ((df>0.25) || (de>0.5));
00637 
00638     if (isBad) ibest=-1;
00639 
00640   }
00641   return ibest;
00642 }
00643