#include <RecoParticleFlow/PFTracking/interface/PFElecTkProducer.h>
Public Member Functions | |
PFElecTkProducer (const edm::ParameterSet &) | |
Constructor. | |
~PFElecTkProducer () | |
Destructor. | |
Private Member Functions | |
virtual void | beginJob (const edm::EventSetup &) |
virtual void | endJob () |
int | FindPfRef (const reco::PFRecTrackCollection &PfRTkColl, reco::GsfTrack, bool) |
bool | otherElId (const reco::GsfTrackCollection &GsfColl, reco::GsfTrack GsfTk) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Produce the PFRecTrack collection. | |
Private Attributes | |
edm::ParameterSet | conf_ |
edm::InputTag | gsfTrackLabel_ |
bool | modemomentum_ |
reco::GsfPFRecTrack | pftrack_ |
edm::InputTag | pfTrackLabel_ |
PFTrackTransformer * | pfTransformer_ |
PFTrackTransformer. | |
bool | trajinev_ |
Trajectory of GSfTracks in the event? |
Definition at line 28 of file PFElecTkProducer.h.
PFElecTkProducer::PFElecTkProducer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Constructor.
Definition at line 33 of file PFElecTkProducer.cc.
References edm::ParameterSet::getParameter(), gsfTrackLabel_, modemomentum_, pfTrackLabel_, and trajinev_.
00033 : 00034 conf_(iConfig), 00035 pfTransformer_(0) 00036 { 00037 LogInfo("PFElecTkProducer")<<"PFElecTkProducer started"; 00038 00039 gsfTrackLabel_ = iConfig.getParameter<InputTag> 00040 ("GsfTrackModuleLabel"); 00041 00042 pfTrackLabel_ = iConfig.getParameter<InputTag> 00043 ("PFRecTrackLabel"); 00044 00045 produces<GsfPFRecTrackCollection>(); 00046 00047 trajinev_ = iConfig.getParameter<bool>("TrajInEvents"); 00048 modemomentum_ = iConfig.getParameter<bool>("ModeMomentum"); 00049 }
PFElecTkProducer::~PFElecTkProducer | ( | ) |
Destructor.
Definition at line 52 of file PFElecTkProducer.cc.
References pfTransformer_.
00053 { 00054 00055 delete pfTransformer_; 00056 }
void PFElecTkProducer::beginJob | ( | const edm::EventSetup & | iSetup | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 291 of file PFElecTkProducer.cc.
References edm::EventSetup::get(), and pfTransformer_.
00292 { 00293 ESHandle<MagneticField> magneticField; 00294 iSetup.get<IdealMagneticFieldRecord>().get(magneticField); 00295 pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0)))); 00296 }
int PFElecTkProducer::FindPfRef | ( | const reco::PFRecTrackCollection & | PfRTkColl, | |
reco::GsfTrack | gsftk, | |||
bool | otherColl | |||
) | [private] |
Definition at line 195 of file PFElecTkProducer.cc.
References reco::TrackBase::eta(), reco::PFRecTrack::KF_ELCAND, muonGeometry::mag(), reco::TrackBase::phi(), reco::Track::seedRef(), funct::sqrt(), and TwoPi.
Referenced by produce().
00197 { 00198 00199 00200 if (&(*gsftk.seedRef())==0) return -1; 00201 00202 reco::PFRecTrackCollection::const_iterator pft=PfRTkColl.begin(); 00203 reco::PFRecTrackCollection::const_iterator pftend=PfRTkColl.end(); 00204 uint i_pf=0; 00205 int ibest=-1; 00206 uint ish_max=0; 00207 float dr_min=1000; 00208 00209 for(;pft!=pftend;++pft){ 00210 uint ish=0; 00211 if ((pft->algoType()==reco::PFRecTrack::KF_ELCAND) || otherColl){ 00212 00213 float dph= fabs(pft->trackRef()->phi()-gsftk.phi()); 00214 if (dph>TMath::TwoPi()) dph-= TMath::TwoPi(); 00215 float det=fabs(pft->trackRef()->eta()-gsftk.eta()); 00216 float dr =sqrt(dph*dph+det*det); 00217 00218 trackingRecHit_iterator hhit= 00219 pft->trackRef()->recHitsBegin(); 00220 trackingRecHit_iterator hhit_end= 00221 pft->trackRef()->recHitsEnd(); 00222 00223 00224 00225 for(;hhit!=hhit_end;++hhit){ 00226 if (!(*hhit)->isValid()) continue; 00227 TrajectorySeed::const_iterator hit= 00228 gsftk.seedRef()->recHits().first; 00229 TrajectorySeed::const_iterator hit_end= 00230 gsftk.seedRef()->recHits().second; 00231 for(;hit!=hit_end;++hit){ 00232 if (!(hit->isValid())) continue; 00233 00234 00235 if((hit->geographicalId()==(*hhit)->geographicalId())&& 00236 (((*hhit)->localPosition()-hit->localPosition()).mag()<0.01)) ish++; 00237 } 00238 00239 } 00240 00241 00242 if ((ish>ish_max)|| 00243 ((ish==ish_max)&&(dr<dr_min))){ 00244 ish_max=ish; 00245 dr_min=dr; 00246 ibest=i_pf; 00247 } 00248 00249 } 00250 00251 i_pf++; 00252 } 00253 if (ibest<0) return -1; 00254 00255 if((ish_max==0) &&(dr_min>0.05))return -1; 00256 if(otherColl && (ish_max==0)) return -1; 00257 return ibest; 00258 }
bool PFElecTkProducer::otherElId | ( | const reco::GsfTrackCollection & | GsfColl, | |
reco::GsfTrack | GsfTk | |||
) | [private] |
Definition at line 261 of file PFElecTkProducer.cc.
References muonGeometry::mag(), reco::TrackBase::numberOfValidHits(), reco::Track::recHitsBegin(), and reco::Track::recHitsEnd().
Referenced by produce().
00262 { 00263 int nhits=GsfTk.numberOfValidHits(); 00264 GsfTrackCollection::const_iterator igs=GsfColl.begin(); 00265 GsfTrackCollection::const_iterator igs_end=GsfColl.end(); 00266 uint shared=0; 00267 for(;igs!=igs_end;++igs){ 00268 uint tmp_sh=0; 00269 trackingRecHit_iterator ghit=igs->recHitsBegin(); 00270 trackingRecHit_iterator ghit_end=igs->recHitsEnd(); 00271 for (;ghit!=ghit_end;++ghit){ 00272 00273 if (!((*ghit)->isValid())) continue; 00274 00275 trackingRecHit_iterator hit=GsfTk.recHitsBegin(); 00276 trackingRecHit_iterator hit_end=GsfTk.recHitsEnd(); 00277 00278 for (;hit!=hit_end;++hit){ 00279 if (!((*hit)->isValid())) continue; 00280 if(((*hit)->geographicalId()==(*ghit)->geographicalId())&& 00281 (((*hit)->localPosition()-(*ghit)->localPosition()).mag()<0.01)) tmp_sh++; 00282 } 00283 } 00284 if (tmp_sh>shared) shared=tmp_sh; 00285 } 00286 return ((float(shared)/float(nhits))<0.5); 00287 }
void PFElecTkProducer::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [private, virtual] |
Produce the PFRecTrack collection.
Implements edm::EDProducer.
Definition at line 65 of file PFElecTkProducer.cc.
References PFTrackTransformer::addPointsAndBrems(), conf_, diff, edm::EventID::event(), FindPfRef(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), reco::PFRecTrack::GSF, gsfTrackLabel_, edm::Event::id(), LogDebug, modemomentum_, otherElId(), pftrack_, pfTrackLabel_, pfTransformer_, edm::Handle< T >::product(), edm::Event::put(), edm::EventID::run(), and trajinev_.
00066 { 00067 LogDebug("PFElecTkProducer")<<"START event: "<<iEvent.id().event() 00068 <<" in run "<<iEvent.id().run(); 00069 00070 //create the empty collections 00071 auto_ptr< GsfPFRecTrackCollection > 00072 gsfPFRecTrackCollection(new GsfPFRecTrackCollection); 00073 00074 00075 //read collections of tracks 00076 Handle<GsfTrackCollection> gsfelectrons; 00077 iEvent.getByLabel(gsfTrackLabel_,gsfelectrons); 00078 00079 //read collections of trajectories 00080 Handle<vector<Trajectory> > TrajectoryCollection; 00081 00082 //read pfrectrack collection 00083 Handle<PFRecTrackCollection> thePfRecTrackCollection; 00084 iEvent.getByLabel(pfTrackLabel_,thePfRecTrackCollection); 00085 const PFRecTrackCollection PfRTkColl = *(thePfRecTrackCollection.product()); 00086 00087 00088 00089 00090 if (trajinev_){ 00091 iEvent.getByLabel(gsfTrackLabel_,TrajectoryCollection); 00092 GsfTrackCollection gsftracks = *(gsfelectrons.product()); 00093 vector<Trajectory> tjvec= *(TrajectoryCollection.product()); 00094 00095 for (uint igsf=0; igsf<gsftracks.size();igsf++) { 00096 00097 GsfTrackRef trackRef(gsfelectrons, igsf); 00098 int kf_ind=FindPfRef(PfRTkColl,gsftracks[igsf],false); 00099 00100 if (kf_ind>=0) { 00101 00102 PFRecTrackRef kf_ref(thePfRecTrackCollection, 00103 kf_ind); 00104 00105 00106 pftrack_=GsfPFRecTrack( gsftracks[igsf].charge(), 00107 reco::PFRecTrack::GSF, 00108 igsf, trackRef, 00109 kf_ref); 00110 } else { 00111 PFRecTrackRef dummyRef; 00112 pftrack_=GsfPFRecTrack( gsftracks[igsf].charge(), 00113 reco::PFRecTrack::GSF, 00114 igsf, trackRef, 00115 dummyRef); 00116 } 00117 00118 bool validgsfbrem = pfTransformer_->addPointsAndBrems(pftrack_, 00119 gsftracks[igsf], 00120 tjvec[igsf], 00121 modemomentum_); 00122 00123 00124 if(validgsfbrem) 00125 gsfPFRecTrackCollection->push_back(pftrack_); 00126 } 00127 //OTHER GSF TRACK COLLECTION 00128 if(conf_.getParameter<bool>("AddGSFTkColl")){ 00129 00130 Handle<GsfElectronCollection> ElecCollection; 00131 iEvent.getByLabel(conf_.getParameter<InputTag >("GsfElectrons"), ElecCollection); 00132 GsfElectronCollection::const_iterator iel = ElecCollection->begin(); 00133 GsfElectronCollection::const_iterator iel_end = ElecCollection->end(); 00134 00135 Handle<GsfTrackCollection> otherGsfColl; 00136 iEvent.getByLabel(conf_.getParameter<InputTag >("GsfTracks"),otherGsfColl); 00137 GsfTrackCollection othergsftracks = *(otherGsfColl.product()); 00138 00139 Handle<vector<Trajectory> > TrajectoryCollection; 00140 iEvent.getByLabel(conf_.getParameter<InputTag >("GsfTracks"),TrajectoryCollection); 00141 vector<Trajectory> newtj= *(TrajectoryCollection.product()); 00142 00143 for(;iel!=iel_end;++iel){ 00144 uint ibest =9999; float diffbest=10000.; 00145 for (uint igsf=0; igsf<othergsftracks.size();igsf++) { 00146 float diff =(iel->gsfTrack()->momentum()-othergsftracks[igsf].momentum()).Mag2(); 00147 if (diff<diffbest){ 00148 ibest=igsf; 00149 diffbest=diff; 00150 } 00151 } 00152 00153 if (ibest==9999 || diffbest>0.00001) continue; 00154 00155 if(otherElId(gsftracks,othergsftracks[ibest])){ 00156 GsfTrackRef trackRef(otherGsfColl, ibest); 00157 00158 int kf_ind=FindPfRef(PfRTkColl,othergsftracks[ibest],true); 00159 00160 if (kf_ind>=0) { 00161 PFRecTrackRef kf_ref(thePfRecTrackCollection, 00162 kf_ind); 00163 pftrack_=GsfPFRecTrack( othergsftracks[ibest].charge(), 00164 reco::PFRecTrack::GSF, 00165 ibest, trackRef, 00166 kf_ref); 00167 } else { 00168 PFRecTrackRef dummyRef; 00169 pftrack_=GsfPFRecTrack( othergsftracks[ibest].charge(), 00170 reco::PFRecTrack::GSF, 00171 ibest, trackRef, 00172 dummyRef); 00173 } 00174 bool validgsfbrem = pfTransformer_->addPointsAndBrems(pftrack_, 00175 othergsftracks[ibest], 00176 newtj[ibest], 00177 modemomentum_); 00178 if(validgsfbrem) 00179 gsfPFRecTrackCollection->push_back(pftrack_); 00180 00181 } 00182 } 00183 } 00184 00185 00186 00187 00188 00189 iEvent.put(gsfPFRecTrackCollection); 00190 }else LogError("PFEleTkProducer")<<"No trajectory in the events"; 00191 00192 }
edm::ParameterSet PFElecTkProducer::conf_ [private] |
bool PFElecTkProducer::modemomentum_ [private] |
edm::InputTag PFElecTkProducer::pfTrackLabel_ [private] |
Definition at line 58 of file PFElecTkProducer.h.
Referenced by beginJob(), produce(), and ~PFElecTkProducer().
bool PFElecTkProducer::trajinev_ [private] |
Trajectory of GSfTracks in the event?
Definition at line 61 of file PFElecTkProducer.h.
Referenced by PFElecTkProducer(), and produce().