Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include "SUSYBSMAnalysis/HSCP/interface/HSCParticleProducer.h"
00023
00024 using namespace susybsm;
00025 using namespace reco;
00026
00027 HSCParticleProducer::HSCParticleProducer(const edm::ParameterSet& iConfig) {
00028 using namespace edm;
00029 using namespace std;
00030
00031
00032 Filter_ = iConfig.getParameter<bool> ("filter");
00033
00034
00035 m_trackTag = iConfig.getParameter<edm::InputTag>("tracks");
00036 m_muonsTag = iConfig.getParameter<edm::InputTag>("muons");
00037 m_trackIsoTag = iConfig.getParameter<edm::InputTag>("tracksIsolation");
00038
00039 useBetaFromTk = iConfig.getParameter<bool> ("useBetaFromTk" );
00040 useBetaFromMuon = iConfig.getParameter<bool> ("useBetaFromMuon");
00041 useBetaFromRpc = iConfig.getParameter<bool> ("useBetaFromRpc" );
00042 useBetaFromEcal = iConfig.getParameter<bool> ("useBetaFromEcal");
00043
00044
00045 minTkP = iConfig.getParameter<double> ("minTkP");
00046 maxTkChi2 = iConfig.getParameter<double> ("maxTkChi2");
00047 minTkHits = iConfig.getParameter<uint32_t>("minTkHits");
00048 minMuP = iConfig.getParameter<double> ("minMuP");
00049 minDR = iConfig.getParameter<double> ("minDR");
00050 maxInvPtDiff = iConfig.getParameter<double> ("maxInvPtDiff");
00051
00052 if(useBetaFromTk )beta_calculator_TK = new BetaCalculatorTK (iConfig);
00053 if(useBetaFromMuon)beta_calculator_MUON = new BetaCalculatorMUON(iConfig);
00054 if(useBetaFromRpc )beta_calculator_RPC = new BetaCalculatorRPC (iConfig);
00055 if(useBetaFromEcal)beta_calculator_ECAL = new BetaCalculatorECAL(iConfig);
00056
00057
00058 std::vector<edm::ParameterSet> SelectionParameters = iConfig.getParameter<std::vector<edm::ParameterSet> >("SelectionParameters");
00059 for(unsigned int i=0;i<SelectionParameters.size();i++){
00060 Selectors.push_back(new CandidateSelector(SelectionParameters[i]) );
00061 }
00062
00063
00064 produces<susybsm::HSCParticleCollection >();
00065 if(useBetaFromEcal)produces<susybsm::HSCPCaloInfoCollection >();
00066
00067 }
00068
00069 HSCParticleProducer::~HSCParticleProducer() {
00070
00071
00072 }
00073
00074
00075
00076
00077
00078
00079 bool
00080 HSCParticleProducer::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00081
00082 using namespace edm;
00083 using namespace reco;
00084 using namespace std;
00085 using namespace susybsm;
00086
00087
00088 edm::Handle<reco::MuonCollection> muonCollectionHandle;
00089 iEvent.getByLabel(m_muonsTag,muonCollectionHandle);
00090
00091
00092 edm::Handle<reco::TrackCollection> trackCollectionHandle;
00093 iEvent.getByLabel(m_trackTag,trackCollectionHandle);
00094
00095
00096 edm::Handle<reco::TrackCollection> trackIsoCollectionHandle;
00097 iEvent.getByLabel(m_trackIsoTag,trackIsoCollectionHandle);
00098
00099
00100
00101 susybsm::HSCParticleCollection* hscp = new susybsm::HSCParticleCollection;
00102 std::auto_ptr<susybsm::HSCParticleCollection> result(hscp);
00103
00104 susybsm::HSCPCaloInfoCollection* caloInfoColl = new susybsm::HSCPCaloInfoCollection;
00105 std::auto_ptr<susybsm::HSCPCaloInfoCollection> caloInfoCollaptr(caloInfoColl);
00106
00107
00108
00109 *hscp = getHSCPSeedCollection(trackCollectionHandle, muonCollectionHandle);
00110
00111
00112 for(susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin(); hscpcandidate != hscp->end(); ++hscpcandidate) {
00113
00114 reco::TrackRef track = hscpcandidate->trackRef();
00115 if(track.isNull())continue;
00116 float dRMin=1000; int found = -1;
00117 for(unsigned int t=0; t<trackIsoCollectionHandle->size();t++) {
00118 reco::TrackRef Isotrack = reco::TrackRef( trackIsoCollectionHandle, t );
00119 if( fabs( (1.0/track->pt())-(1.0/Isotrack->pt())) > maxInvPtDiff) continue;
00120 float dR = deltaR(track->momentum(), Isotrack->momentum());
00121 if(dR <= minDR && dR < dRMin){ dRMin=dR; found = t;}
00122 }
00123 if(found>=0)hscpcandidate->setTrackIso(reco::TrackRef( trackIsoCollectionHandle, found ));
00124 }
00125
00126
00127 if(useBetaFromTk){
00128 for(susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin(); hscpcandidate != hscp->end(); ++hscpcandidate) {
00129 beta_calculator_TK->addInfoToCandidate(*hscpcandidate, iEvent,iSetup);
00130 }}
00131
00132
00133 if(useBetaFromMuon){
00134 for(susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin(); hscpcandidate != hscp->end(); ++hscpcandidate) {
00135 beta_calculator_MUON->addInfoToCandidate(*hscpcandidate, iEvent,iSetup);
00136 }}
00137
00138
00139 if(useBetaFromRpc){
00140 for(susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin(); hscpcandidate != hscp->end(); ++hscpcandidate) {
00141 beta_calculator_RPC->addInfoToCandidate(*hscpcandidate, iEvent, iSetup);
00142 }}
00143
00144
00145
00146
00147
00148 if(useBetaFromEcal){
00149 int Index=0;
00150 caloInfoColl->resize(hscp->size());
00151 for(susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin(); hscpcandidate != hscp->end(); ++hscpcandidate, Index++) {
00152 beta_calculator_ECAL->addInfoToCandidate(*hscpcandidate,trackCollectionHandle,iEvent,iSetup, (*caloInfoColl)[Index]);
00153 }}
00154
00155
00156 for(int i=0;i<(int)hscp->size();i++) {
00157 susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin() + i;
00158 bool decision = false;
00159 for(unsigned int s=0;s<Selectors.size();s++){decision |= Selectors[s]->isSelected(*hscpcandidate);}
00160 if(!decision){
00161 hscp->erase(hscpcandidate);
00162 if(useBetaFromEcal)caloInfoColl->erase(caloInfoColl->begin() + i);
00163 i--;
00164 }
00165 }
00166 bool filterResult = !Filter_ || (Filter_ && hscp->size()>=1);
00167
00168
00169
00170
00171
00172 if(useBetaFromEcal){
00173 edm::OrphanHandle<susybsm::HSCPCaloInfoCollection> caloInfoHandle= iEvent.put(caloInfoCollaptr);
00174
00175 for(int i=0;i<(int)hscp->size();i++) {
00176 susybsm::HSCParticleCollection::iterator hscpcandidate = hscp->begin() + i;
00177 hscpcandidate->setCaloInfo(HSCPCaloInfoRef(caloInfoHandle,i));
00178 }
00179 }
00180
00181
00182
00183
00184
00185 edm::OrphanHandle<susybsm::HSCParticleCollection> putHandle = iEvent.put(result);
00186
00187
00188
00189
00190
00191
00192
00193 return filterResult;
00194 }
00195
00196
00197 void
00198 HSCParticleProducer::beginJob() {
00199 }
00200
00201
00202 void
00203 HSCParticleProducer::endJob() {
00204 }
00205
00206 std::vector<HSCParticle> HSCParticleProducer::getHSCPSeedCollection(edm::Handle<reco::TrackCollection>& trackCollectionHandle, edm::Handle<reco::MuonCollection>& muonCollectionHandle)
00207 {
00208 std::vector<HSCParticle> HSCPCollection;
00209
00210
00211 std::vector<reco::TrackRef> tracks;
00212 for(unsigned int i=0; i<trackCollectionHandle->size(); i++){
00213 TrackRef track = reco::TrackRef( trackCollectionHandle, i );
00214 if(track->p()<minTkP || (track->chi2()/track->ndof())>maxTkChi2 || track->found()<minTkHits)continue;
00215 tracks.push_back( track );
00216 }
00217
00218
00219 for(unsigned int m=0; m<muonCollectionHandle->size(); m++){
00220 reco::MuonRef muon = reco::MuonRef( muonCollectionHandle, m );
00221 if(muon->p()<minMuP )continue;
00222 TrackRef innertrack = muon->innerTrack();
00223 if(innertrack.isNull())continue;
00224
00225
00226
00227 float dRMin=1000; int found = -1;
00228 for(unsigned int t=0; t<tracks.size();t++) {
00229 reco::TrackRef track = tracks[t];
00230 if( fabs( (1.0/innertrack->pt())-(1.0/track->pt())) > maxInvPtDiff) continue;
00231 float dR = deltaR(innertrack->momentum(), track->momentum());
00232 if(dR <= minDR && dR < dRMin){ dRMin=dR; found = t;}
00233 }
00234
00235 HSCParticle candidate;
00236 candidate.setMuon(muon);
00237 if(found>=0){
00238
00239 candidate.setTrack(tracks[found]);
00240 tracks.erase(tracks.begin()+found);
00241 }
00242 HSCPCollection.push_back(candidate);
00243 }
00244
00245
00246 for(unsigned int m=0; m<muonCollectionHandle->size(); m++){
00247 reco::MuonRef muon = reco::MuonRef( muonCollectionHandle, m );
00248 if(muon->p()<minMuP)continue;
00249 TrackRef innertrack = muon->innerTrack();
00250 if(innertrack.isNonnull())continue;
00251
00252
00253 float dRMin=1000; int found = -1;
00254 for(unsigned int t=0; t<tracks.size();t++) {
00255 reco::TrackRef track = tracks[t];
00256 if( fabs( (1.0/muon->pt())-(1.0/track->pt())) > maxInvPtDiff) continue;
00257 float dR = deltaR(muon->momentum(), track->momentum());
00258 if(dR <= minDR && dR < dRMin){ dRMin=dR; found = t;}
00259 }
00260
00261 HSCParticle candidate;
00262 candidate.setMuon(muon);
00263 if(found>=0){
00264
00265 candidate.setTrack(tracks[found]);
00266 tracks.erase(tracks.begin()+found);
00267 }
00268 HSCPCollection.push_back(candidate);
00269 }
00270
00271
00272 for(unsigned int i=0; i<tracks.size(); i++){
00273 HSCParticle candidate;
00274 candidate.setTrack(tracks[i]);
00275 HSCPCollection.push_back(candidate);
00276 }
00277
00278 return HSCPCollection;
00279 }
00280
00281
00282 DEFINE_FWK_MODULE(HSCParticleProducer);
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301