-
Notifications
You must be signed in to change notification settings - Fork 4.2k
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Merge pull request #8135 from wmtford/clusterRefinerTagMCmerged
Test module that tags merged clusters per MC truth
- Loading branch information
Showing
6 changed files
with
179 additions
and
1 deletion.
There are no files selected for viewing
19 changes: 19 additions & 0 deletions
19
RecoLocalTracker/SiStripClusterizer/python/test/ClusterRefinerTagMCmerged_cfi.py
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,19 @@ | ||
import FWCore.ParameterSet.Config as cms | ||
|
||
siStripClusters = cms.EDProducer("ClusterRefinerTagMCmerged", | ||
UntaggedClusterProducer = cms.InputTag('siStripClustersUntagged'), | ||
ClusterRefiner = cms.PSet( | ||
# For TrackerHitAssociator | ||
ROUList = cms.vstring('g4SimHitsTrackerHitsTIBLowTof', | ||
'g4SimHitsTrackerHitsTIBHighTof', | ||
'g4SimHitsTrackerHitsTIDLowTof', | ||
'g4SimHitsTrackerHitsTIDHighTof', | ||
'g4SimHitsTrackerHitsTOBLowTof', | ||
'g4SimHitsTrackerHitsTOBHighTof', | ||
'g4SimHitsTrackerHitsTECLowTof', | ||
'g4SimHitsTrackerHitsTECHighTof'), | ||
associateRecoTracks = cms.bool(True), # True to save some time if no PU | ||
associatePixel = cms.bool(False), | ||
associateStrip = cms.bool(True) | ||
) | ||
) |
26 changes: 26 additions & 0 deletions
26
RecoLocalTracker/SiStripClusterizer/python/test/tagMCmergedCustomize_cff.py
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,26 @@ | ||
# | ||
# With this customization the ClusterizerRefinerTagMCmerged module will be substituted | ||
# for the standard clusterizer. If a cluster is matched to more than one simTrack | ||
# its "merged" bit will be set, so that SiStripCluster::isMerged() will return true. | ||
# | ||
# If pileup is present, add the following line so that only in-time simTracks will | ||
# be counted, and make sure that process.mix is on the path. | ||
# process.siStripClusters.ClusterRefiner.associateRecoTracks = cms.bool(False) | ||
# | ||
|
||
import FWCore.ParameterSet.Config as cms | ||
|
||
def tagMCmerged(process): | ||
|
||
process.siStripClustersUntagged = process.siStripClusters.clone() | ||
stripClusIndex = process.striptrackerlocalreco.index(process.siStripClusters) | ||
process.striptrackerlocalreco.remove(process.siStripClusters) | ||
del process.siStripClusters | ||
process.load('RecoLocalTracker.SiStripClusterizer.test.ClusterRefinerTagMCmerged_cfi') | ||
process.siStripClustersTagged = cms.Sequence(process.siStripClustersUntagged*process.siStripClusters) | ||
process.striptrackerlocalreco.insert(stripClusIndex,process.siStripClustersTagged) | ||
|
||
# Override the chargePerCM cut in stripCPE and use cluster::isMerged() instead. | ||
process.StripCPEfromTrackAngleESProducer.parameters.maxChgOneMIP = cms.double(-6000.) | ||
|
||
return(process) |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -1,4 +1,5 @@ | ||
<library name="RecoLocalTrackerSiStripClusterizerTestPlugins" file="*.cc"> | ||
<use name="RecoLocalTracker/SiStripClusterizer"/> | ||
<use name="SimTracker/TrackerHitAssociation"/> | ||
<flags EDM_PLUGIN="1"/> | ||
</library> |
90 changes: 90 additions & 0 deletions
90
RecoLocalTracker/SiStripClusterizer/test/ClusterRefinerTagMCmerged.cc
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,90 @@ | ||
#include "RecoLocalTracker/SiStripClusterizer/test/ClusterRefinerTagMCmerged.h" | ||
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h" | ||
#include "FWCore/Framework/interface/Event.h" | ||
|
||
ClusterRefinerTagMCmerged:: | ||
ClusterRefinerTagMCmerged(const edm::ParameterSet& conf) | ||
: inputTag( conf.getParameter<edm::InputTag>("UntaggedClusterProducer") ), | ||
confClusterRefiner_(conf.getParameter<edm::ParameterSet>("ClusterRefiner")) { | ||
produces< edmNew::DetSetVector<SiStripCluster> > (); | ||
inputToken = consumes< edmNew::DetSetVector<SiStripCluster> >(inputTag); | ||
} | ||
|
||
void ClusterRefinerTagMCmerged:: | ||
produce(edm::Event& event, const edm::EventSetup& es) { | ||
|
||
std::auto_ptr< edmNew::DetSetVector<SiStripCluster> > output(new edmNew::DetSetVector<SiStripCluster>()); | ||
output->reserve(10000,4*10000); | ||
|
||
associator_.reset(new TrackerHitAssociator(event, confClusterRefiner_)); | ||
edm::Handle< edmNew::DetSetVector<SiStripCluster> > input; | ||
|
||
if ( findInput(inputToken, input, event) ) refineCluster(input, output); | ||
else edm::LogError("Input Not Found") << "[ClusterRefinerTagMCmerged::produce] ";// << inputTag; | ||
|
||
LogDebug("Output") << output->dataSize() << " clusters from " | ||
<< output->size() << " modules"; | ||
output->shrink_to_fit(); | ||
event.put(output); | ||
} | ||
|
||
void ClusterRefinerTagMCmerged:: | ||
refineCluster(const edm::Handle< edmNew::DetSetVector<SiStripCluster> >& input, | ||
std::auto_ptr< edmNew::DetSetVector<SiStripCluster> >& output) { | ||
|
||
for (edmNew::DetSetVector<SiStripCluster>::const_iterator det=input->begin(); det!=input->end(); det++) { | ||
// DetSetVector filler to receive the clusters we produce | ||
edmNew::DetSetVector<SiStripCluster>::FastFiller outFill(*output, det->id()); | ||
uint32_t detId = det->id(); | ||
int ntk = 0; | ||
int NtkAll = 0; | ||
for (edmNew::DetSet<SiStripCluster>::iterator clust = det->begin(); clust != det->end(); clust++) { | ||
std::vector<uint8_t> amp = clust->amplitudes(); | ||
SiStripCluster* newCluster = new SiStripCluster(clust->firstStrip(), amp.begin(), amp.end()); | ||
if (associator_ != 0) { | ||
std::vector<SimHitIdpr> simtrackid; | ||
bool useAssociateHit = !confClusterRefiner_.getParameter<bool>("associateRecoTracks"); | ||
if (useAssociateHit) { | ||
std::vector<PSimHit> simhit; | ||
associator_->associateCluster(clust, DetId(detId), simtrackid, simhit); | ||
NtkAll = simtrackid.size(); | ||
ntk = 0; | ||
if (simtrackid.size() > 1) { | ||
for (auto const& it : simtrackid) { | ||
int NintimeHit = 0; | ||
for (auto const& ih : simhit) { | ||
// std::cout << " hit(tk, evt) trk(tk, evt) bunch (" | ||
// << ih.trackId() << ", " << ih.eventId().rawId() << ") (" | ||
// << it.first << ", " << it.second.rawId() << ") " | ||
// << ih.eventId().bunchCrossing() | ||
// << std::endl; | ||
if (ih.trackId() == it.first && ih.eventId() == it.second && ih.eventId().bunchCrossing() == 0) ++NintimeHit; | ||
} | ||
if (NintimeHit > 0) ++ntk; | ||
} | ||
} | ||
} else { | ||
associator_->associateSimpleRecHitCluster(clust, DetId(detId), simtrackid); | ||
ntk = NtkAll = simtrackid.size(); | ||
} | ||
if (ntk > 1) { | ||
newCluster->setMerged(true); | ||
} else { | ||
newCluster->setMerged(false); | ||
} | ||
} | ||
outFill.push_back(*newCluster); | ||
// std::cout << "t_m_w " << " " << NtkAll << " " << newCluster->isMerged() | ||
// << " " << clust->amplitudes().size() | ||
// << std::endl; | ||
} // traverse clusters | ||
} // traverse sensors | ||
} | ||
|
||
template<class T> | ||
inline | ||
bool ClusterRefinerTagMCmerged:: | ||
findInput(const edm::EDGetTokenT<T>& tag, edm::Handle<T>& handle, const edm::Event& e) { | ||
e.getByToken( tag, handle); | ||
return handle.isValid(); | ||
} |
41 changes: 41 additions & 0 deletions
41
RecoLocalTracker/SiStripClusterizer/test/ClusterRefinerTagMCmerged.h
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,41 @@ | ||
#ifndef RecoLocalTracker_ClusterRefinerTagMCmerged_h | ||
#define RecoLocalTracker_ClusterRefinerTagMCmerged_h | ||
|
||
// | ||
// Use MC truth to identify merged clusters, i.e., those associated with more than one | ||
// (in-time) SimTrack. | ||
// | ||
// Author: Bill Ford (wtford) 6 March 2015 | ||
// | ||
|
||
#include "FWCore/Framework/interface/Frameworkfwd.h" | ||
#include "FWCore/Framework/interface/stream/EDProducer.h" | ||
#include "FWCore/ParameterSet/interface/ParameterSet.h" | ||
#include "FWCore/Utilities/interface/InputTag.h" | ||
#include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h" | ||
|
||
#include <memory> | ||
|
||
class ClusterRefinerTagMCmerged : public edm::stream::EDProducer<> { | ||
|
||
public: | ||
|
||
explicit ClusterRefinerTagMCmerged(const edm::ParameterSet& conf); | ||
virtual void produce(edm::Event&, const edm::EventSetup&); | ||
|
||
private: | ||
|
||
template<class T> bool findInput(const edm::EDGetTokenT<T>&, edm::Handle<T>&, const edm::Event&); | ||
virtual void refineCluster(const edm::Handle< edmNew::DetSetVector<SiStripCluster> >& input, | ||
std::auto_ptr< edmNew::DetSetVector<SiStripCluster> >& output); | ||
|
||
const edm::InputTag inputTag; | ||
typedef edm::EDGetTokenT< edmNew::DetSetVector<SiStripCluster> > token_t; | ||
token_t inputToken; | ||
edm::ParameterSet confClusterRefiner_; | ||
|
||
std::shared_ptr<TrackerHitAssociator> associator_; | ||
|
||
}; | ||
|
||
#endif |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters