Go to the documentation of this file.00001
00002
00003
00004 #include "PhysicsTools/HepMCCandAlgos/interface/FlavorHistoryProducer.h"
00005
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "DataFormats/Math/interface/deltaR.h"
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034 using namespace std;
00035 using namespace reco;
00036 using namespace edm;
00037
00038 FlavorHistoryProducer::FlavorHistoryProducer( const ParameterSet & p ) :
00039 src_( p.getParameter<InputTag>( "src" ) ),
00040 matchedSrc_( p.getParameter<InputTag>( "matchedSrc") ),
00041 matchDR_ ( p.getParameter<double> ("matchDR") ),
00042 pdgIdToSelect_( p.getParameter<int> ("pdgIdToSelect") ),
00043 ptMinParticle_( p.getParameter<double>( "ptMinParticle") ),
00044 ptMinShower_( p.getParameter<double>( "ptMinShower") ),
00045 etaMaxParticle_( p.getParameter<double>( "etaMaxParticle" )),
00046 etaMaxShower_( p.getParameter<double>( "etaMaxShower" )),
00047 flavorHistoryName_( p.getParameter<string>("flavorHistoryName") ),
00048 verbose_( p.getUntrackedParameter<bool>( "verbose" ) )
00049 {
00050 produces<FlavorHistoryEvent >(flavorHistoryName_);
00051 }
00052
00053 FlavorHistoryProducer::~FlavorHistoryProducer() {
00054 }
00055
00056 void FlavorHistoryProducer::produce( Event& evt, const EventSetup& )
00057 {
00058 if ( verbose_ ) cout << "---------- Hello from FlavorHistoryProducer! -----" << endl;
00059
00060
00061 Handle<CandidateView > particlesViewH;
00062 evt.getByLabel( src_, particlesViewH );
00063
00064 Handle<CandidateView> matchedView;
00065 evt.getByLabel( matchedSrc_, matchedView );
00066
00067
00068 vector<const Candidate* > particles;
00069 for( CandidateView::const_iterator p = particlesViewH->begin(); p != particlesViewH->end(); ++p ) {
00070 particles.push_back(&*p);
00071 }
00072
00073
00074 vector<int> partonIndices;
00075
00076 vector<int> progenitorIndices;
00077
00078 vector<int> sisterIndices;
00079
00080 vector<FlavorHistory::FLAVOR_T> flavorSources;
00081
00082
00083 auto_ptr<FlavorHistoryEvent > flavorHistoryEvent ( new FlavorHistoryEvent () ) ;
00084
00085
00086
00087
00088 vector<const Candidate* >::size_type j;
00089 vector<const Candidate* >::size_type j_max = particles.size();
00090 for( j=0; j<j_max; ++j ) {
00091
00092
00093 const Candidate *p = particles[j];
00094
00095 vector<Candidate const *>::size_type partonIndex=j;
00096 vector<Candidate const *>::size_type progenitorIndex=j_max;
00097 bool foundProgenitor = false;
00098 FlavorHistory::FLAVOR_T flavorSource=FlavorHistory::FLAVOR_NULL;
00099
00100
00101 int idabs = std::abs( (p)->pdgId() );
00102 int nDa = (p)->numberOfDaughters();
00103
00104
00105
00106 if ( idabs <= FlavorHistory::tQuarkId && p->status() == 2 &&
00107 idabs == pdgIdToSelect_ ) {
00108
00109 if ( nDa > 0 ) {
00110
00111 if((p)->pt() > ptMinShower_ && fabs((p)->eta())<etaMaxShower_) {
00112
00113 if(verbose_) cout << "--------------------------" << endl;
00114 if(verbose_) cout << "Processing particle " << j << ", status = " << p->status() << ", pdg id = " << p->pdgId() << endl;
00115
00116
00117
00118 vector<Candidate const *> allParents;
00119 getAncestors( *p, allParents );
00120
00121 if(verbose_) {
00122 cout << "Parents : " << endl;
00123 vector<Candidate const *>::const_iterator iprint = allParents.begin(),
00124 iprintend = allParents.end();
00125 for ( ; iprint != iprintend; ++iprint )
00126 cout << " status = " << (*iprint)->status() << ", pdg id = " << (*iprint)->pdgId() << ", pt = " << (*iprint)->pt() << endl;
00127 }
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149 vector<Candidate const *>::size_type a_size = allParents.size();
00150
00151
00152
00153
00154 for( vector<Candidate const *>::size_type i=0 ; i < a_size && !foundProgenitor; ++i ) {
00155 const Candidate * aParent=allParents[i];
00156 vector<Candidate const *>::const_iterator found = find(particles.begin(),particles.end(),aParent);
00157
00158
00159
00160 progenitorIndex = found - particles.begin();
00161
00162 int aParentId = std::abs(aParent->pdgId());
00163
00164
00165
00166
00167 if( progenitorIndex > 5 ) {
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177 if( aParent->numberOfMothers() == 2 &&
00178 aParent->pdgId() == p->pdgId() && aParent->status() == 3 ) {
00179 if(verbose_) cout << "Matrix Element progenitor" << endl;
00180 if ( flavorSource == FlavorHistory::FLAVOR_NULL ) flavorSource = FlavorHistory::FLAVOR_ME;
00181 foundProgenitor = true;
00182 }
00183
00184 else if( (aParentId>pdgIdToSelect_ && aParentId<FlavorHistory::gluonId) ||
00185 aParentId > FlavorHistory::gluonId ) {
00186 if(verbose_) cout << "Flavor decay progenitor" << endl;
00187 if ( flavorSource == FlavorHistory::FLAVOR_NULL ) flavorSource = FlavorHistory::FLAVOR_DECAY;
00188 foundProgenitor = true;
00189 }
00190 }
00191 }
00192
00193
00194
00195 if ( !foundProgenitor ) {
00196 if ( flavorSource == FlavorHistory::FLAVOR_NULL ) flavorSource = FlavorHistory::FLAVOR_GS;
00197
00198 for( vector<Candidate const *>::size_type i=0 ; i < a_size && !foundProgenitor; ++i ) {
00199 const Candidate * aParent=allParents[i];
00200 vector<Candidate const *>::const_iterator found = find(particles.begin(),particles.end(),aParent);
00201
00202 progenitorIndex = found - particles.begin();
00203
00204 if ( aParent->numberOfMothers() == 1 && aParent->status() == 3 ) {
00205 foundProgenitor = true;
00206 }
00207 }
00208 }
00209
00210 }
00211 }
00212
00213
00214
00215
00216 if ( !foundProgenitor ) progenitorIndex = j_max;
00217
00218
00219
00220 partonIndices.push_back( partonIndex );
00221 progenitorIndices.push_back( progenitorIndex );
00222 flavorSources.push_back(flavorSource);
00223 sisterIndices.push_back( -1 );
00224
00225 }
00226 }
00227
00228
00229
00230
00231
00232 if ( verbose_ ) cout << "Making sisters" << endl;
00233
00234 if ( partonIndices.size() == progenitorIndices.size() && partonIndices.size() > 0 ) {
00235
00236 for ( unsigned int ii = 0; ii < partonIndices.size(); ++ii ) {
00237
00238 const Candidate * candi = particles[partonIndices[ii]];
00239 if ( verbose_ ) cout << " Sister candidate particle 1: " << ii << ", pdgid = " << candi->pdgId() << ", status = " << candi->status() << endl;
00240
00241
00242
00243 for ( unsigned int jj = 0; jj < partonIndices.size(); ++jj ) {
00244 if ( ii != jj ) {
00245 const Candidate * candj = particles[partonIndices[jj]];
00246 if ( verbose_ ) cout << " Sister candidate particle 2: " << jj << ", pdgid = " << candj->pdgId() << ", status = " << candj->status() << endl;
00247
00248
00249 if ( candi->pdgId() == -1 * candj->pdgId() && candi->status() == candj->status() ) {
00250 sisterIndices[ii] = partonIndices[jj];
00251 if ( verbose_ ) cout << "Parton " << partonIndices[ii] << " has sister " << sisterIndices[ii] << endl;
00252 }
00253 }
00254 }
00255
00256
00257 if ( sisterIndices[ii] < 0 ) {
00258 if ( verbose_ ) cout << "No sister. Classified as flavor excitation" << endl;
00259 flavorSources[ii] = FlavorHistory::FLAVOR_EXC;
00260 }
00261
00262 if ( verbose_ ) cout << "Getting matched jet" << endl;
00263
00264 CandidateView::const_iterator matched = getClosestJet( matchedView, *candi );
00265 CandidateView::const_iterator matchedBegin = matchedView->begin();
00266 CandidateView::const_iterator matchedEnd = matchedView->end();
00267
00268
00269 if ( verbose_ ) cout << "Getting sister jet" << endl;
00270
00271 CandidateView::const_iterator sister =
00272 ( sisterIndices[ii] >= 0 && static_cast<unsigned int>(sisterIndices[ii]) < particles.size() ) ?
00273 getClosestJet( matchedView, *particles[sisterIndices[ii]] ) :
00274 matchedEnd ;
00275
00276 if ( verbose_ ) cout << "Making matched shallow clones : ";
00277 ShallowClonePtrCandidate matchedCand ;
00278 if ( matched != matchedEnd ) {
00279 if ( verbose_ ) cout << " found" << endl;
00280 matchedCand = ShallowClonePtrCandidate( CandidatePtr(matchedView, matched - matchedBegin ) );
00281 } else {
00282 if ( verbose_ ) cout << " NOT found" << endl;
00283 }
00284
00285 if ( verbose_ ) cout << "Making sister shallow clones : ";
00286 ShallowClonePtrCandidate sisterCand;
00287 if ( sister != matchedEnd ) {
00288 if ( verbose_ ) cout << " found" << endl;
00289 sisterCand = ShallowClonePtrCandidate( CandidatePtr(matchedView, sister - matchedBegin ) );
00290 } else {
00291 if ( verbose_ ) cout << " NOT found" << endl;
00292 }
00293
00294 if ( verbose_ ) cout << "Making history object" << endl;
00295
00296 FlavorHistory history (flavorSources[ii],
00297 particlesViewH,
00298 partonIndices[ii], progenitorIndices[ii], sisterIndices[ii],
00299 matchedCand,
00300 sisterCand );
00301 if ( verbose_ ) cout << "Adding flavor history : " << history << endl;
00302 flavorHistoryEvent->push_back( history );
00303 }
00304 }
00305
00306
00307 if ( flavorHistoryEvent->size() > 0 ) {
00308
00309
00310 if ( verbose_ ) cout << "Caching flavor history event" << endl;
00311 flavorHistoryEvent->cache();
00312
00313 if ( verbose_ ) {
00314 cout << "Outputting pdg id = " << pdgIdToSelect_ << " with nelements = " << flavorHistoryEvent->size() << endl;
00315 vector<FlavorHistory>::const_iterator i = flavorHistoryEvent->begin(),
00316 iend = flavorHistoryEvent->end();
00317 for ( ; i !=iend; ++i ) {
00318 cout << *i << endl;
00319 }
00320 }
00321 }
00322
00323
00324 evt.put( flavorHistoryEvent, flavorHistoryName_ );
00325 }
00326
00327
00328
00329 void FlavorHistoryProducer::getAncestors(const Candidate &c,
00330 vector<Candidate const *> & moms )
00331 {
00332
00333 if( c.numberOfMothers() == 1 ) {
00334 const Candidate * dau = &c;
00335 const Candidate * mom = c.mother();
00336 while ( dau->numberOfMothers() != 0) {
00337 moms.push_back( dau );
00338 dau = mom ;
00339 mom = dau->mother();
00340 }
00341 }
00342 }
00343
00344
00345 CandidateView::const_iterator
00346 FlavorHistoryProducer::getClosestJet( Handle<CandidateView> const & pJets,
00347 reco::Candidate const & parton ) const
00348 {
00349 double dr = matchDR_;
00350 CandidateView::const_iterator j = pJets->begin(),
00351 jend = pJets->end();
00352 CandidateView::const_iterator bestJet = pJets->end();
00353 for ( ; j != jend; ++j ) {
00354 double dri = deltaR( parton.p4(), j->p4() );
00355 if ( dri < dr ) {
00356 dr = dri;
00357 bestJet = j;
00358 }
00359 }
00360 return bestJet;
00361 }
00362
00363 #include "FWCore/Framework/interface/MakerMacros.h"
00364
00365 DEFINE_FWK_MODULE( FlavorHistoryProducer );