Skip to content

Commit

Permalink
Merge pull request #8989 from mark-grimes/addInitialVertices75X
Browse files Browse the repository at this point in the history
Add option to save initial vertex collection (forward port of #8975)
  • Loading branch information
cmsbuild committed May 8, 2015
2 parents a7cfc3e + 430a8ce commit 9b5bea3
Show file tree
Hide file tree
Showing 3 changed files with 53 additions and 3 deletions.
Expand Up @@ -4,6 +4,7 @@
accumulatorType = cms.string('TrackingTruthAccumulator'),
createUnmergedCollection = cms.bool(True),
createMergedBremsstrahlung = cms.bool(True),
createInitialVertexCollection = cms.bool(False),
alwaysAddAncestors = cms.bool(True),
maximumPreviousBunchCrossing = cms.uint32(9999),
maximumSubsequentBunchCrossing = cms.uint32(9999),
Expand Down
45 changes: 42 additions & 3 deletions SimGeneral/TrackingAnalysis/plugins/TrackingTruthAccumulator.cc
Expand Up @@ -15,17 +15,20 @@
* @date circa Oct/2012 to Feb/2013
*
* Changelog:
* 05/May/2015 Mark Grimes - Added functionality to add a collection of just the initial vertices
* for FastTimer studies.
*
* 17/Jul/2014 Dominik Nowatschin (dominik.nowatschin@cern.ch) - added SimVertex and a ref to
* HepMC::Genvertex to TrackingVertex in TrackingParticleFactory::createTrackingVertex; handle to
* edm::HepMCProduct is created directly in TrackingTruthAccumulator::accumulate and not in
* accumulateEvent as edm::Event and PileUpEventPrincipal have different getByLabel() functions
*
* 07/Feb/2013 Mark Grimes - Reorganised and added a bit more documentation. Still not enough
* though.
* 12/Mar/2012 (branch NewTrackingParticle only) Mark Grimes - Updated TrackingParticle creation
* to fit in with Subir Sarkar's re-implementation of TrackingParticle.
* 12/Mar/2012 Mark Grimes - Updated TrackingParticle creation to fit in with Subir Sarkar's
* re-implementation of TrackingParticle.
*/
#include "SimGeneral/TrackingAnalysis/interface/TrackingTruthAccumulator.h"
#include "SimGeneral/TrackingAnalysis/plugins/TrackingTruthAccumulator.h"

#include "SimGeneral/MixingModule/interface/DigiAccumulatorMixModFactory.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
Expand Down Expand Up @@ -224,6 +227,7 @@ TrackingTruthAccumulator::TrackingTruthAccumulator( const edm::ParameterSet & co
maximumSubsequentBunchCrossing_( config.getParameter<unsigned int>("maximumSubsequentBunchCrossing") ),
createUnmergedCollection_( config.getParameter<bool>("createUnmergedCollection") ),
createMergedCollection_(config.getParameter<bool>("createMergedBremsstrahlung") ),
createInitialVertexCollection_(config.getParameter<bool>("createInitialVertexCollection") ),
addAncestors_( config.getParameter<bool>("alwaysAddAncestors") ),
removeDeadModules_( config.getParameter<bool>("removeDeadModules") ),
simTrackLabel_( config.getParameter<edm::InputTag>("simTrackCollection") ),
Expand Down Expand Up @@ -287,6 +291,11 @@ TrackingTruthAccumulator::TrackingTruthAccumulator( const edm::ParameterSet & co
mixMod.produces<TrackingVertexCollection>("MergedTrackTruth");
}

if( createInitialVertexCollection_ )
{
mixMod.produces<TrackingVertexCollection>("InitialVertices");
}

iC.consumes<std::vector<SimTrack> >(simTrackLabel_);
iC.consumes<std::vector<SimVertex> >(simVertexLabel_);
iC.consumes<std::vector<reco::GenParticle> >(genParticleLabel_);
Expand Down Expand Up @@ -326,6 +335,11 @@ void TrackingTruthAccumulator::initializeEvent( edm::Event const& event, edm::Ev
mergedOutput_.refTrackingParticles=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingParticleCollection>("MergedTrackTruth");
mergedOutput_.refTrackingVertexes=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingVertexCollection>("MergedTrackTruth");
}

if( createInitialVertexCollection_ )
{
pInitialVertices_.reset( new TrackingVertexCollection );
}
}

/// create handle to edm::HepMCProduct here because event.getByLabel with edm::HepMCProduct only works for edm::Event
Expand Down Expand Up @@ -378,6 +392,12 @@ void TrackingTruthAccumulator::finalizeEvent( edm::Event& event, edm::EventSetup
event.put( mergedOutput_.pTrackingVertices, "MergedTrackTruth" );
}

if( createInitialVertexCollection_ )
{
edm::LogInfo("TrackingTruthAccumulator") << "Adding " << pInitialVertices_->size() << " initial TrackingVertexs to the event.";

event.put( pInitialVertices_, "InitialVertices" );
}
}

template<class T> void TrackingTruthAccumulator::accumulateEvent( const T& event, const edm::EventSetup& setup, const edm::Handle< edm::HepMCProduct >& hepMCproduct)
Expand Down Expand Up @@ -470,6 +490,25 @@ template<class T> void TrackingTruthAccumulator::accumulateEvent( const T& event
// specifies adding ancestors, the function is called recursively to do that.
::addTrack( pDecayTrack, pSelector, pUnmergedCollectionWrapper.get(), pMergedCollectionWrapper.get(), objectFactory, addAncestors_, tTopo );
}

// If configured to create a collection of initial vertices, add them from this bunch
// crossing. No selection is applied on this collection, but it also has no links to
// the TrackingParticle decay products.
// There are a lot of "initial vertices", I'm not entirely sure what they all are
// (nuclear interactions in the detector maybe?), but the one for the main event is
// the one with vertexId==0.
if( createInitialVertexCollection_ )
{
// Pretty sure the one with vertexId==0 is always the first one, but doesn't hurt to check
for( const auto& pRootVertex : decayChain.rootVertices )
{
const SimVertex& vertex=hSimVertices->at(decayChain.rootVertices[0]->simVertexIndex);
if( vertex.vertexId()!=0 ) continue;

pInitialVertices_->push_back( objectFactory.createTrackingVertex(pRootVertex) );
break;
}
}
}

template<class T> void TrackingTruthAccumulator::fillSimHits( std::vector<const PSimHit*>& returnValue, const T& event, const edm::EventSetup& setup )
Expand Down
Expand Up @@ -46,6 +46,13 @@ class PSimHit;
* <tr><td> createMergedBremsstrahlung </td><td> bool </td><td> Whether to create the TrackingParticle collection with bremsstrahlung merged. At
* least one of createUnmergedCollection or createMergedBremsstrahlung should be true
* otherwise nothing will be produced. </td></tr>
* <tr><td> createInitialVertexCollection </td><td> bool </td><td> Whether to create a collection of just the initial vertices. You can usually get this
* information from one of the other collections (merged or unmerged bremsstrahlung), but
* for this collection no selection is applied. Hence you will always have all of the initial
* vertices regardless of how tightly you select TrackingParticles with the "select" parameter.
* Note that <b>the collection will have no links to the products of these vertices<b>. If you
* want to know what came off these vertices you will have to look in one of the other
* collections. The name of the collection will be "InitialVertices".</td></tr>
* <tr><td> alwaysAddAncestors </td><td> bool </td><td> If a sim track passes selection and is turned into a TrackingParticle, all of it's
* parents will also be created even if they fail the selection. This was the default
* behaviour for the old TrackingParticleProducer. </td></tr>
Expand Down Expand Up @@ -98,6 +105,8 @@ class TrackingTruthAccumulator : public DigiAccumulatorMixMod
/// If bremsstrahlung merging, whether to also add the unmerged collection to the event or not.
const bool createUnmergedCollection_;
const bool createMergedCollection_;
/// Whether or not to create a separate collection for just the initial interaction vertices
const bool createInitialVertexCollection_;
/// Whether or not to add the full parentage of any TrackingParticle that is inserted in the collection.
const bool addAncestors_;

Expand Down Expand Up @@ -138,6 +147,7 @@ class TrackingTruthAccumulator : public DigiAccumulatorMixMod
private:
OutputCollections unmergedOutput_;
OutputCollections mergedOutput_;
std::auto_ptr<TrackingVertexCollection> pInitialVertices_;
};

#endif // end of "#ifndef TrackingAnalysis_TrackingTruthAccumulator_h"

0 comments on commit 9b5bea3

Please sign in to comment.