CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes
MultiTrackValidator Class Reference

#include <MultiTrackValidator.h>

Inheritance diagram for MultiTrackValidator:
DQMEDAnalyzer edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > > edm::stream::EDAnalyzerBase edm::EDConsumerBase MultiTrackValidatorGenPs

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 Method called once per event. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 Method called to book the DQM histograms. More...
 
 MultiTrackValidator (const edm::ParameterSet &pset)
 Constructor. More...
 
virtual ~MultiTrackValidator ()
 Destructor. More...
 
- Public Member Functions inherited from DQMEDAnalyzer
void beginRun (edm::Run const &, edm::EventSetup const &) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer (void)
 
void endLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
void endRunSummary (edm::Run const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
uint32_t streamId () const
 
- Public Member Functions inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
 EDAnalyzer ()=default
 
- Public Member Functions inherited from edm::stream::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Protected Attributes

std::vector< edm::InputTagassociators
 
edm::EDGetTokenT< reco::BeamSpotbsSrc
 
const bool calculateDrSingleCollection_
 
const bool dodEdxPlots_
 
const bool doMVAPlots_
 
const bool doPlotsOnlyForTruePV_
 
const bool doPVAssociationPlots_
 
const bool doRecoTrackPlots_
 
std::vector< bool > doResolutionPlots_
 
const bool doSeedPlots_
 
const bool doSimPlots_
 
const bool doSimTrackPlots_
 
const bool doSummaryPlots_
 
std::unique_ptr< MTVHistoProducerAlgoForTrackerhistoProducerAlgo_
 
const bool ignoremissingtkcollection_
 
std::vector< edm::InputTaglabel
 
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > label_pileupinfo
 
edm::EDGetTokenT< TrackingParticleCollectionlabel_tp_effic
 
edm::EDGetTokenT< TrackingParticleRefVectorlabel_tp_effic_refvector
 
edm::EDGetTokenT< TrackingParticleCollectionlabel_tp_fake
 
edm::EDGetTokenT< TrackingParticleRefVectorlabel_tp_fake_refvector
 
edm::EDGetTokenT< TrackingVertexCollectionlabel_tv
 
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > labelToken
 
std::vector< edm::EDGetTokenT< edm::View< TrajectorySeed > > > labelTokenSeed
 
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx1Tag
 
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx2Tag
 
std::string parametersDefiner
 
const bool parametersDefinerIsCosmic_
 
std::vector< edm::EDGetTokenT< std::vector< PSimHit > > > simHitTokens_
 
const bool useAssociators_
 

Private Types

using MVACollection = std::vector< float >
 
using QualityMaskCollection = std::vector< unsigned char >
 

Private Member Functions

const reco::Vertex::PointgetRecoPVPosition (const edm::Event &event, const edm::Handle< TrackingVertexCollection > &htv) const
 
const TrackingVertex::LorentzVectorgetSimPVPosition (const edm::Handle< TrackingVertexCollection > &htv) const
 
size_t tpDR (const TrackingParticleRefVector &tPCeff, const std::vector< size_t > &selected_tPCeff, DynArray< float > &dR_tPCeff) const
 
void tpParametersAndSelection (const TrackingParticleRefVector &tPCeff, const ParametersDefinerForTP &parametersDefinerTP, const edm::Event &event, const edm::EventSetup &setup, const reco::BeamSpot &bs, std::vector< std::tuple< TrackingParticle::Vector, TrackingParticle::Point > > &momVert_tPCeff, std::vector< size_t > &selected_tPCeff) const
 
void trackDR (const edm::View< reco::Track > &trackCollection, const edm::View< reco::Track > &trackCollectionDr, DynArray< float > &dR_trk) const
 

Private Attributes

edm::EDGetTokenT< SimHitTPAssociationProducer::SimHitTPAssociationList_simHitTpMapTag
 
std::vector< edm::EDGetTokenT< reco::RecoToSimCollection > > associatormapRtSs
 
std::vector< edm::EDGetTokenT< reco::SimToRecoCollection > > associatormapStRs
 
std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associatorTokens
 
CosmicTrackingParticleSelector cosmictpSelector
 
std::string dirName_
 
TrackingParticleSelector dRtpSelector
 
std::unique_ptr< RecoTrackSelectorBasedRTrackSelector
 
std::vector< MonitorElement * > h_assoc2_coll
 
std::vector< MonitorElement * > h_assoc_coll
 
std::vector< MonitorElement * > h_looper_coll
 
std::vector< MonitorElement * > h_pileup_coll
 
std::vector< MonitorElement * > h_reco_coll
 
std::vector< MonitorElement * > h_simul_coll
 
edm::EDGetTokenT< edm::View< reco::Track > > labelTokenForDrCalculation
 
std::vector< std::vector< std::tuple< edm::EDGetTokenT< MVACollection >, edm::EDGetTokenT< QualityMaskCollection > > > > mvaQualityCollectionTokens_
 
edm::EDGetTokenT< edm::View< reco::Vertex > > recoVertexToken_
 
const double simPVMaxZ_
 
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNLayersToken_
 
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNPixelLayersToken_
 
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNStripStereoLayersToken_
 
TrackingParticleSelector tpSelector
 
bool useGsf
 
edm::EDGetTokenT< reco::VertexToTrackingVertexAssociatorvertexAssociatorToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDAnalyzerBase
typedef EDAnalyzerAdaptorBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static std::shared_ptr< dqmDetails::NoCacheglobalBeginLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *)
 
static std::shared_ptr< dqmDetails::NoCacheglobalBeginRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *)
 
static void globalEndLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *, dqmDetails::NoCache *)
 
static void globalEndRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *, dqmDetails::NoCache *)
 
- Static Public Member Functions inherited from edm::stream::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Class that prodecs histrograms to validate Track Reconstruction performances

Author
cerati

Definition at line 32 of file MultiTrackValidator.h.

Member Typedef Documentation

using MultiTrackValidator::MVACollection = std::vector<float>
private

Definition at line 113 of file MultiTrackValidator.h.

using MultiTrackValidator::QualityMaskCollection = std::vector<unsigned char>
private

Definition at line 114 of file MultiTrackValidator.h.

Constructor & Destructor Documentation

MultiTrackValidator::MultiTrackValidator ( const edm::ParameterSet pset)

Constructor.

Definition at line 52 of file MultiTrackValidator.cc.

References _simHitTpMapTag, associatormapRtSs, associatormapStRs, associators, associatorTokens, AlignmentProducer_cff::beamSpotTag, begin, bsSrc, calculateDrSingleCollection_, edm::EDConsumerBase::consumes(), edm::EDConsumerBase::consumesCollector(), cosmictpSelector, dirName_, dodEdxPlots_, doMVAPlots_, doPlotsOnlyForTruePV_, doPVAssociationPlots_, doResolutionPlots_, MultiTrackValidator_cfi::doResolutionPlotsForLabels, doSeedPlots_, dRtpSelector, dRTrackSelector, relativeConstraints::empty, end, Exception, spr::find(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), histoProducerAlgo_, checklumidiff::l, edm::InputTag::label(), label, label_pileupinfo, label_tp_effic, label_tp_effic_refvector, label_tp_fake, label_tp_fake_refvector, label_tv, tablePrinter::labels, edm::EDConsumerBase::labelsForToken(), labelToken, labelTokenForDrCalculation, m_dEdx1Tag, m_dEdx2Tag, MTVHistoProducerAlgoForTracker::makeRecoTrackSelectorFromTPSelectorParameters(), edm::ProductLabels::module, mvaQualityCollectionTokens_, recoVertexToken_, simHitTokens_, TrackRefitter_38T_cff::src, AlCaHLTBitMon_QueryRunRegistry::string, GlobalPosition_Frontier_DevDB_cff::tag, tpNLayersToken_, tpNPixelLayersToken_, tpNStripStereoLayersToken_, tpSelector, useAssociators_, useGsf, edm::vector_transform(), and vertexAssociatorToken_.

52  :
53  associators(pset.getUntrackedParameter< std::vector<edm::InputTag> >("associators")),
54  label(pset.getParameter< std::vector<edm::InputTag> >("label")),
55  parametersDefiner(pset.getParameter<std::string>("parametersDefiner")),
56  parametersDefinerIsCosmic_(parametersDefiner == "CosmicParametersDefinerForTP"),
57  ignoremissingtkcollection_(pset.getUntrackedParameter<bool>("ignoremissingtrackcollection",false)),
58  useAssociators_(pset.getParameter< bool >("UseAssociators")),
59  calculateDrSingleCollection_(pset.getUntrackedParameter<bool>("calculateDrSingleCollection")),
60  doPlotsOnlyForTruePV_(pset.getUntrackedParameter<bool>("doPlotsOnlyForTruePV")),
61  doSummaryPlots_(pset.getUntrackedParameter<bool>("doSummaryPlots")),
62  doSimPlots_(pset.getUntrackedParameter<bool>("doSimPlots")),
63  doSimTrackPlots_(pset.getUntrackedParameter<bool>("doSimTrackPlots")),
64  doRecoTrackPlots_(pset.getUntrackedParameter<bool>("doRecoTrackPlots")),
65  dodEdxPlots_(pset.getUntrackedParameter<bool>("dodEdxPlots")),
66  doPVAssociationPlots_(pset.getUntrackedParameter<bool>("doPVAssociationPlots")),
67  doSeedPlots_(pset.getUntrackedParameter<bool>("doSeedPlots")),
68  doMVAPlots_(pset.getUntrackedParameter<bool>("doMVAPlots")),
69  simPVMaxZ_(pset.getUntrackedParameter<double>("simPVMaxZ"))
70 {
71  const edm::InputTag& label_tp_effic_tag = pset.getParameter< edm::InputTag >("label_tp_effic");
72  const edm::InputTag& label_tp_fake_tag = pset.getParameter< edm::InputTag >("label_tp_fake");
73 
74  if(pset.getParameter<bool>("label_tp_effic_refvector")) {
75  label_tp_effic_refvector = consumes<TrackingParticleRefVector>(label_tp_effic_tag);
76  }
77  else {
78  label_tp_effic = consumes<TrackingParticleCollection>(label_tp_effic_tag);
79  }
80  if(pset.getParameter<bool>("label_tp_fake_refvector")) {
81  label_tp_fake_refvector = consumes<TrackingParticleRefVector>(label_tp_fake_tag);
82  }
83  else {
84  label_tp_fake = consumes<TrackingParticleCollection>(label_tp_fake_tag);
85  }
86  label_pileupinfo = consumes<std::vector<PileupSummaryInfo> >(pset.getParameter< edm::InputTag >("label_pileupinfo"));
87  for(const auto& tag: pset.getParameter<std::vector<edm::InputTag>>("sim")) {
88  simHitTokens_.push_back(consumes<std::vector<PSimHit>>(tag));
89  }
90 
91  std::vector<edm::InputTag> doResolutionPlotsForLabels = pset.getParameter<std::vector<edm::InputTag> >("doResolutionPlotsForLabels");
92  doResolutionPlots_.reserve(label.size());
93  for (auto& itag : label) {
94  labelToken.push_back(consumes<edm::View<reco::Track> >(itag));
95  const bool doResol = doResolutionPlotsForLabels.empty() || (std::find(cbegin(doResolutionPlotsForLabels), cend(doResolutionPlotsForLabels), itag) != cend(doResolutionPlotsForLabels));
96  doResolutionPlots_.push_back(doResol);
97  }
98  { // check for duplicates
99  auto labelTmp = edm::vector_transform(label, [&](const edm::InputTag& tag) { return tag.label(); });
100  std::sort(begin(labelTmp), end(labelTmp));
102  const std::string* prev = &empty;
103  for(const std::string& l: labelTmp) {
104  if(l == *prev) {
105  throw cms::Exception("Configuration") << "Duplicate InputTag in labels: " << l;
106  }
107  prev = &l;
108  }
109  }
110 
112  bsSrc = consumes<reco::BeamSpot>(beamSpotTag);
113 
114  ParameterSet psetForHistoProducerAlgo = pset.getParameter<ParameterSet>("histoProducerAlgoBlock");
115  histoProducerAlgo_ = std::make_unique<MTVHistoProducerAlgoForTracker>(psetForHistoProducerAlgo, beamSpotTag, doSeedPlots_, consumesCollector());
116 
117  dirName_ = pset.getParameter<std::string>("dirName");
118 
119  tpNLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_nlayers"));
120  tpNPixelLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_npixellayers"));
121  tpNStripStereoLayersToken_ = consumes<edm::ValueMap<unsigned int> >(pset.getParameter<edm::InputTag>("label_tp_nstripstereolayers"));
122 
123  if(dodEdxPlots_) {
124  m_dEdx1Tag = consumes<edm::ValueMap<reco::DeDxData> >(pset.getParameter< edm::InputTag >("dEdx1Tag"));
125  m_dEdx2Tag = consumes<edm::ValueMap<reco::DeDxData> >(pset.getParameter< edm::InputTag >("dEdx2Tag"));
126  }
127 
128  label_tv = consumes<TrackingVertexCollection>(pset.getParameter< edm::InputTag >("label_tv"));
130  recoVertexToken_ = consumes<edm::View<reco::Vertex> >(pset.getUntrackedParameter<edm::InputTag>("label_vertex"));
131  vertexAssociatorToken_ = consumes<reco::VertexToTrackingVertexAssociator>(pset.getUntrackedParameter<edm::InputTag>("vertexAssociator"));
132  }
133 
134  if(doMVAPlots_) {
136  auto mvaPSet = pset.getUntrackedParameter<edm::ParameterSet>("mvaLabels");
137  for(size_t iIter=0; iIter<labelToken.size(); ++iIter) {
139  labelsForToken(labelToken[iIter], labels);
140  if(mvaPSet.exists(labels.module)) {
141  mvaQualityCollectionTokens_[iIter] = edm::vector_transform(mvaPSet.getUntrackedParameter<std::vector<std::string> >(labels.module),
142  [&](const std::string& tag) {
143  return std::make_tuple(consumes<MVACollection>(edm::InputTag(tag, "MVAValues")),
144  consumes<QualityMaskCollection>(edm::InputTag(tag, "QualityMasks")));
145  });
146  }
147  }
148  }
149 
150  tpSelector = TrackingParticleSelector(pset.getParameter<double>("ptMinTP"),
151  pset.getParameter<double>("ptMaxTP"),
152  pset.getParameter<double>("minRapidityTP"),
153  pset.getParameter<double>("maxRapidityTP"),
154  pset.getParameter<double>("tipTP"),
155  pset.getParameter<double>("lipTP"),
156  pset.getParameter<int>("minHitTP"),
157  pset.getParameter<bool>("signalOnlyTP"),
158  pset.getParameter<bool>("intimeOnlyTP"),
159  pset.getParameter<bool>("chargedOnlyTP"),
160  pset.getParameter<bool>("stableOnlyTP"),
161  pset.getParameter<std::vector<int> >("pdgIdTP"));
162 
164  pset.getParameter<double>("minRapidityTP"),
165  pset.getParameter<double>("maxRapidityTP"),
166  pset.getParameter<double>("tipTP"),
167  pset.getParameter<double>("lipTP"),
168  pset.getParameter<int>("minHitTP"),
169  pset.getParameter<bool>("chargedOnlyTP"),
170  pset.getParameter<std::vector<int> >("pdgIdTP"));
171 
172 
173  ParameterSet psetVsPhi = psetForHistoProducerAlgo.getParameter<ParameterSet>("TpSelectorForEfficiencyVsPhi");
174  dRtpSelector = TrackingParticleSelector(psetVsPhi.getParameter<double>("ptMin"),
175  psetVsPhi.getParameter<double>("ptMax"),
176  psetVsPhi.getParameter<double>("minRapidity"),
177  psetVsPhi.getParameter<double>("maxRapidity"),
178  psetVsPhi.getParameter<double>("tip"),
179  psetVsPhi.getParameter<double>("lip"),
180  psetVsPhi.getParameter<int>("minHit"),
181  psetVsPhi.getParameter<bool>("signalOnly"),
182  psetVsPhi.getParameter<bool>("intimeOnly"),
183  psetVsPhi.getParameter<bool>("chargedOnly"),
184  psetVsPhi.getParameter<bool>("stableOnly"),
185  psetVsPhi.getParameter<std::vector<int> >("pdgId"));
186 
188 
189  useGsf = pset.getParameter<bool>("useGsf");
190 
191  _simHitTpMapTag = mayConsume<SimHitTPAssociationProducer::SimHitTPAssociationList>(pset.getParameter<edm::InputTag>("simHitTpMapTag"));
192 
194  labelTokenForDrCalculation = consumes<edm::View<reco::Track> >(pset.getParameter<edm::InputTag>("trackCollectionForDrCalculation"));
195  }
196 
197  if(useAssociators_) {
198  for (auto const& src: associators) {
199  associatorTokens.push_back(consumes<reco::TrackToTrackingParticleAssociator>(src));
200  }
201  } else {
202  for (auto const& src: associators) {
203  associatormapStRs.push_back(consumes<reco::SimToRecoCollection>(src));
204  associatormapRtSs.push_back(consumes<reco::RecoToSimCollection>(src));
205  }
206  }
207 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< SimHitTPAssociationProducer::SimHitTPAssociationList > _simHitTpMapTag
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx2Tag
std::vector< edm::EDGetTokenT< reco::SimToRecoCollection > > associatormapStRs
CosmicTrackingParticleSelector cosmictpSelector
std::vector< std::vector< std::tuple< edm::EDGetTokenT< MVACollection >, edm::EDGetTokenT< QualityMaskCollection > > > > mvaQualityCollectionTokens_
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNStripStereoLayersToken_
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNLayersToken_
TrackingParticleSelector dRtpSelector
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associatorTokens
const bool doPlotsOnlyForTruePV_
edm::EDGetTokenT< reco::BeamSpot > bsSrc
static std::unique_ptr< RecoTrackSelectorBase > makeRecoTrackSelectorFromTPSelectorParameters(const edm::ParameterSet &pset, const edm::InputTag &beamSpotTag, edm::ConsumesCollector &iC)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNPixelLayersToken_
std::vector< edm::InputTag > associators
SingleObjectSelector< TrackingParticleCollection,::TrackingParticleSelector > TrackingParticleSelector
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > labelToken
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx1Tag
std::vector< edm::EDGetTokenT< reco::RecoToSimCollection > > associatormapRtSs
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
TrackingParticleSelector tpSelector
std::vector< bool > doResolutionPlots_
edm::EDGetTokenT< edm::View< reco::Track > > labelTokenForDrCalculation
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > label_pileupinfo
std::vector< edm::InputTag > label
#define end
Definition: vmac.h:37
char const * module
Definition: ProductLabels.h:5
ObjectSelector< CosmicTrackingParticleSelector > CosmicTrackingParticleSelector
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_fake_refvector
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
edm::EDGetTokenT< TrackingVertexCollection > label_tv
edm::EDGetTokenT< reco::VertexToTrackingVertexAssociator > vertexAssociatorToken_
std::unique_ptr< RecoTrackSelectorBase > dRTrackSelector
const bool ignoremissingtkcollection_
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
std::string const & label() const
Definition: InputTag.h:36
const bool doPVAssociationPlots_
edm::EDGetTokenT< edm::View< reco::Vertex > > recoVertexToken_
std::vector< edm::EDGetTokenT< std::vector< PSimHit > > > simHitTokens_
const bool calculateDrSingleCollection_
#define begin
Definition: vmac.h:30
edm::EDGetTokenT< TrackingParticleCollection > label_tp_fake
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_effic_refvector
doResolutionPlotsForLabels
do resolution plots only for these labels (or all if empty)
edm::EDGetTokenT< TrackingParticleCollection > label_tp_effic
const bool parametersDefinerIsCosmic_
MultiTrackValidator::~MultiTrackValidator ( )
virtual

Destructor.

Definition at line 210 of file MultiTrackValidator.cc.

210 {}

Member Function Documentation

void MultiTrackValidator::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
override

Method called once per event.

Definition at line 482 of file MultiTrackValidator.cc.

References _simHitTpMapTag, funct::abs(), reco::TrackToTrackingParticleAssociator::associateRecoToSim(), reco::TrackToTrackingParticleAssociator::associateSimToReco(), associationMapFilterValues(), associatormapRtSs, associatormapStRs, associators, associatorTokens, bsSrc, EncodedEventId::bunchCrossing(), calculateDrSingleCollection_, ParametersDefinerForTP::clone(), cosmictpSelector, declareDynArray, dodEdxPlots_, doMVAPlots_, doPlotsOnlyForTruePV_, doPVAssociationPlots_, doRecoTrackPlots_, doResolutionPlots_, doSeedPlots_, doSimTrackPlots_, doSummaryPlots_, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, dRtpSelector, dRTrackSelector, TrackingParticleIP::dxy(), TrackingParticleIP::dz(), relativeConstraints::empty, edm::AssociationMap< Tag >::end(), EncodedEventId::event(), TrackingParticle::eventId(), Exception, edm::AssociationMap< Tag >::find(), plotBeamSpotDB::first, edm::EventSetup::get(), edm::Event::getByToken(), PileupSummaryInfo::getPU_NumInteractions(), getRecoPVPosition(), getSimPVPosition(), h_assoc2_coll, h_assoc_coll, h_looper_coll, h_pileup_coll, h_reco_coll, h_simul_coll, reco::TrackBase::highPurity, histoProducerAlgo_, mps_fire::i, ignoremissingtkcollection_, CosmicTrackingParticleSelector::initEvent(), edm::EDGetTokenT< T >::isUninitialized(), label, label_pileupinfo, label_tp_effic, label_tp_effic_refvector, label_tp_fake, label_tp_fake_refvector, label_tv, labelToken, labelTokenForDrCalculation, LogDebug, LogTrace, reco::TrackBase::loose, m_dEdx1Tag, m_dEdx2Tag, SiStripPI::max, TrackingParticle::momentum(), eostools::move(), DetachedQuadStep_cff::mva, mvaQualityCollectionTokens_, TrackingParticle::numberOfTrackerHits(), parametersDefiner, parametersDefinerIsCosmic_, reco::BeamSpot::position(), edm::Handle< T >::product(), edm::RefToBaseVector< T >::push_back(), edm::RefVector< C, T, F >::push_back(), edm::View< T >::refAt(), simPVMaxZ_, edm::RefVector< C, T, F >::size(), edm::View< T >::size(), findQualityFiles::size, mathSSE::sqrt(), tpDR(), tpNLayersToken_, tpNPixelLayersToken_, tpNStripStereoLayersToken_, tpParametersAndSelection(), HiIsolationCommonParameters_cff::track, findElectronsInSiStrips_cfi::trackCollection, TrackValidationHeavyIons_cff::trackCollectionForDrCalculation, trackDR(), trackFromSeedFitFailed(), useAssociators_, TrackingParticle::vertex(), and w.

482  {
483  using namespace reco;
484 
485  LogDebug("TrackValidator") << "\n====================================================" << "\n"
486  << "Analyzing new event" << "\n"
487  << "====================================================\n" << "\n";
488 
489 
490  edm::ESHandle<ParametersDefinerForTP> parametersDefinerTPHandle;
491  setup.get<TrackAssociatorRecord>().get(parametersDefiner,parametersDefinerTPHandle);
492  //Since we modify the object, we must clone it
493  auto parametersDefinerTP = parametersDefinerTPHandle->clone();
494 
496  setup.get<TrackerTopologyRcd>().get(httopo);
497  const TrackerTopology& ttopo = *httopo;
498 
499  // FIXME: we really need to move to edm::View for reading the
500  // TrackingParticles... Unfortunately it has non-trivial
501  // consequences on the associator/association interfaces etc.
502  TrackingParticleRefVector tmpTPeff;
503  TrackingParticleRefVector tmpTPfake;
504  const TrackingParticleRefVector *tmpTPeffPtr = nullptr;
505  const TrackingParticleRefVector *tmpTPfakePtr = nullptr;
506 
508  edm::Handle<TrackingParticleRefVector> TPCollectionHeffRefVector;
509 
510  const bool tp_effic_refvector = label_tp_effic.isUninitialized();
511  if(!tp_effic_refvector) {
512  event.getByToken(label_tp_effic, TPCollectionHeff);
513  for(size_t i=0, size=TPCollectionHeff->size(); i<size; ++i) {
514  tmpTPeff.push_back(TrackingParticleRef(TPCollectionHeff, i));
515  }
516  tmpTPeffPtr = &tmpTPeff;
517  }
518  else {
519  event.getByToken(label_tp_effic_refvector, TPCollectionHeffRefVector);
520  tmpTPeffPtr = TPCollectionHeffRefVector.product();
521  }
523  edm::Handle<TrackingParticleCollection> TPCollectionHfake ;
524  event.getByToken(label_tp_fake,TPCollectionHfake);
525  for(size_t i=0, size=TPCollectionHfake->size(); i<size; ++i) {
526  tmpTPfake.push_back(TrackingParticleRef(TPCollectionHfake, i));
527  }
528  tmpTPfakePtr = &tmpTPfake;
529  }
530  else {
531  edm::Handle<TrackingParticleRefVector> TPCollectionHfakeRefVector;
532  event.getByToken(label_tp_fake_refvector, TPCollectionHfakeRefVector);
533  tmpTPfakePtr = TPCollectionHfakeRefVector.product();
534  }
535 
536  TrackingParticleRefVector const & tPCeff = *tmpTPeffPtr;
537  TrackingParticleRefVector const & tPCfake = *tmpTPfakePtr;
538 
539  ensureEffIsSubsetOfFake(tPCeff, tPCfake);
540 
543  //warning: make sure the TP collection used in the map is the same used in the MTV!
544  event.getByToken(_simHitTpMapTag,simHitsTPAssoc);
545  parametersDefinerTP->initEvent(simHitsTPAssoc);
546  cosmictpSelector.initEvent(simHitsTPAssoc);
547  }
548  dRTrackSelector->init(event, setup);
549  histoProducerAlgo_->init(event, setup);
550 
551  // Find the sim PV and tak its position
553  event.getByToken(label_tv, htv);
554  const TrackingVertex::LorentzVector *theSimPVPosition = getSimPVPosition(htv);
555  if(simPVMaxZ_ >= 0) {
556  if(!theSimPVPosition) return;
557  if(std::abs(theSimPVPosition->z()) > simPVMaxZ_) return;
558  }
559 
560  // Check, when necessary, if reco PV matches to sim PV
561  const reco::Vertex::Point *thePVposition = nullptr;
563  thePVposition = getRecoPVPosition(event, htv);
564  if(doPlotsOnlyForTruePV_ && !thePVposition)
565  return;
566 
567  // Rest of the code assumes that if thePVposition is non-null, the
568  // PV-association histograms get filled. In above, the "nullness"
569  // is used to deliver the information if the reco PV is matched to
570  // the sim PV.
572  thePVposition = nullptr;
573  }
574 
575  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
576  event.getByToken(bsSrc,recoBeamSpotHandle);
577  reco::BeamSpot const & bs = *recoBeamSpotHandle;
578 
580  event.getByToken(label_pileupinfo,puinfoH);
581  PileupSummaryInfo puinfo;
582 
583  for (unsigned int puinfo_ite=0;puinfo_ite<(*puinfoH).size();++puinfo_ite){
584  if ((*puinfoH)[puinfo_ite].getBunchCrossing()==0){
585  puinfo=(*puinfoH)[puinfo_ite];
586  break;
587  }
588  }
589 
590  // Number of 3D layers for TPs
592  event.getByToken(tpNLayersToken_, tpNLayersH);
593  const auto& nLayers_tPCeff = *tpNLayersH;
594 
595  event.getByToken(tpNPixelLayersToken_, tpNLayersH);
596  const auto& nPixelLayers_tPCeff = *tpNLayersH;
597 
598  event.getByToken(tpNStripStereoLayersToken_, tpNLayersH);
599  const auto& nStripMonoAndStereoLayers_tPCeff = *tpNLayersH;
600 
601  // Precalculate TP selection (for efficiency), and momentum and vertex wrt PCA
602  //
603  // TODO: ParametersDefinerForTP ESProduct needs to be changed to
604  // EDProduct because of consumes.
605  //
606  // In principle, we could just precalculate the momentum and vertex
607  // wrt PCA for all TPs for once and put that to the event. To avoid
608  // repetitive calculations those should be calculated only once for
609  // each TP. That would imply that we should access TPs via Refs
610  // (i.e. View) in here, since, in general, the eff and fake TP
611  // collections can be different (and at least HI seems to use that
612  // feature). This would further imply that the
613  // RecoToSimCollection/SimToRecoCollection should be changed to use
614  // View<TP> instead of vector<TP>, and migrate everything.
615  //
616  // Or we could take only one input TP collection, and do another
617  // TP-selection to obtain the "fake" collection like we already do
618  // for "efficiency" TPs.
619  std::vector<size_t> selected_tPCeff;
620  std::vector<std::tuple<TrackingParticle::Vector, TrackingParticle::Point>> momVert_tPCeff;
621  tpParametersAndSelection(tPCeff, *parametersDefinerTP, event, setup, bs, momVert_tPCeff, selected_tPCeff);
622 
623  //calculate dR for TPs
624  declareDynArray(float, tPCeff.size(), dR_tPCeff);
625  size_t n_selTP_dr = tpDR(tPCeff, selected_tPCeff, dR_tPCeff);
626 
629  event.getByToken(labelTokenForDrCalculation, trackCollectionForDrCalculation);
630  }
631 
632  // dE/dx
633  // at some point this could be generalized, with a vector of tags and a corresponding vector of Handles
634  // I'm writing the interface such to take vectors of ValueMaps
635  std::vector<const edm::ValueMap<reco::DeDxData> *> v_dEdx;
636  if(dodEdxPlots_) {
639  event.getByToken(m_dEdx1Tag, dEdx1Handle);
640  event.getByToken(m_dEdx2Tag, dEdx2Handle);
641  v_dEdx.push_back(dEdx1Handle.product());
642  v_dEdx.push_back(dEdx2Handle.product());
643  }
644 
645  std::vector<const MVACollection *> mvaCollections;
646  std::vector<const QualityMaskCollection *> qualityMaskCollections;
647  std::vector<float> mvaValues;
648 
649  int w=0; //counter counting the number of sets of histograms
650  for (unsigned int ww=0;ww<associators.size();ww++){
651  // run value filtering of recoToSim map already here as it depends only on the association, not track collection
652  reco::SimToRecoCollection const * simRecCollPFull=nullptr;
653  reco::RecoToSimCollection const * recSimCollP=nullptr;
654  reco::RecoToSimCollection recSimCollL;
655  if(!useAssociators_) {
656  Handle<reco::SimToRecoCollection > simtorecoCollectionH;
657  event.getByToken(associatormapStRs[ww], simtorecoCollectionH);
658  simRecCollPFull = simtorecoCollectionH.product();
659 
660  Handle<reco::RecoToSimCollection > recotosimCollectionH;
661  event.getByToken(associatormapRtSs[ww],recotosimCollectionH);
662  recSimCollP = recotosimCollectionH.product();
663 
664  // We need to filter the associations of the fake-TrackingParticle
665  // collection only from RecoToSim collection, otherwise the
666  // RecoToSim histograms get false entries
667  recSimCollL = associationMapFilterValues(*recSimCollP, tPCfake);
668  recSimCollP = &recSimCollL;
669  }
670 
671  for (unsigned int www=0;www<label.size();www++, w++){ // need to increment w here, since there are many continues in the loop body
672  //
673  //get collections from the event
674  //
675  edm::Handle<View<Track> > trackCollectionHandle;
676  if(!event.getByToken(labelToken[www], trackCollectionHandle)&&ignoremissingtkcollection_)continue;
677  const edm::View<Track>& trackCollection = *trackCollectionHandle;
678 
679  reco::SimToRecoCollection const * simRecCollP=nullptr;
680  reco::SimToRecoCollection simRecCollL;
681 
682  //associate tracks
683  LogTrace("TrackValidator") << "Analyzing "
684  << label[www] << " with "
685  << associators[ww] <<"\n";
686  if(useAssociators_){
688  event.getByToken(associatorTokens[ww], theAssociator);
689 
690  // The associator interfaces really need to be fixed...
692  for(edm::View<Track>::size_type i=0; i<trackCollection.size(); ++i) {
693  trackRefs.push_back(trackCollection.refAt(i));
694  }
695 
696 
697  LogTrace("TrackValidator") << "Calling associateRecoToSim method" << "\n";
698  recSimCollL = std::move(theAssociator->associateRecoToSim(trackRefs, tPCfake));
699  recSimCollP = &recSimCollL;
700  LogTrace("TrackValidator") << "Calling associateSimToReco method" << "\n";
701  // It is necessary to do the association wrt. fake TPs,
702  // because this SimToReco association is used also for
703  // duplicates. Since the set of efficiency TPs are required to
704  // be a subset of the set of fake TPs, for efficiency
705  // histograms it doesn't matter if the association contains
706  // associations of TPs not in the set of efficiency TPs.
707  simRecCollL = std::move(theAssociator->associateSimToReco(trackRefs, tPCfake));
708  simRecCollP = &simRecCollL;
709  }
710  else{
711  // We need to filter the associations of the current track
712  // collection only from SimToReco collection, otherwise the
713  // SimToReco histograms get false entries. The filtering must
714  // be done separately for each track collection.
715  simRecCollL = associationMapFilterValues(*simRecCollPFull, trackCollection);
716  simRecCollP = &simRecCollL;
717  }
718 
719  reco::RecoToSimCollection const & recSimColl = *recSimCollP;
720  reco::SimToRecoCollection const & simRecColl = *simRecCollP;
721 
722  // read MVA collections
726  for(const auto& tokenTpl: mvaQualityCollectionTokens_[www]) {
727  event.getByToken(std::get<0>(tokenTpl), hmva);
728  event.getByToken(std::get<1>(tokenTpl), hqual);
729 
730  mvaCollections.push_back(hmva.product());
731  qualityMaskCollections.push_back(hqual.product());
732  if(mvaCollections.back()->size() != trackCollection.size()) {
733  throw cms::Exception("Configuration") << "Inconsistency in track collection and MVA sizes. Track collection " << www << " has " << trackCollection.size() << " tracks, whereas the MVA " << (mvaCollections.size()-1) << " for it has " << mvaCollections.back()->size() << " entries. Double-check your configuration.";
734  }
735  if(qualityMaskCollections.back()->size() != trackCollection.size()) {
736  throw cms::Exception("Configuration") << "Inconsistency in track collection and quality mask sizes. Track collection " << www << " has " << trackCollection.size() << " tracks, whereas the quality mask " << (qualityMaskCollections.size()-1) << " for it has " << qualityMaskCollections.back()->size() << " entries. Double-check your configuration.";
737  }
738  }
739  }
740 
741  // ########################################################
742  // fill simulation histograms (LOOP OVER TRACKINGPARTICLES)
743  // ########################################################
744 
745  //compute number of tracks per eta interval
746  //
747  LogTrace("TrackValidator") << "\n# of TrackingParticles: " << tPCeff.size() << "\n";
748  int ats(0); //This counter counts the number of simTracks that are "associated" to recoTracks
749  int st(0); //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
750 
751  //loop over already-selected TPs for tracking efficiency
752  for(size_t i=0; i<selected_tPCeff.size(); ++i) {
753  size_t iTP = selected_tPCeff[i];
754  const TrackingParticleRef& tpr = tPCeff[iTP];
755  const TrackingParticle& tp = *tpr;
756 
757  auto const& momVert = momVert_tPCeff[i];
758  TrackingParticle::Vector momentumTP;
759  TrackingParticle::Point vertexTP;
760 
761  double dxySim(0);
762  double dzSim(0);
763  double dxyPVSim = 0;
764  double dzPVSim = 0;
765  double dR=dR_tPCeff[iTP];
766 
767  //---------- THIS PART HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
768  //If the TrackingParticle is collison like, get the momentum and vertex at production state
770  {
771  momentumTP = tp.momentum();
772  vertexTP = tp.vertex();
773  //Calcualte the impact parameters w.r.t. PCA
774  const TrackingParticle::Vector& momentum = std::get<TrackingParticle::Vector>(momVert);
775  const TrackingParticle::Point& vertex = std::get<TrackingParticle::Point>(momVert);
776  dxySim = TrackingParticleIP::dxy(vertex, momentum, bs.position());
777  dzSim = TrackingParticleIP::dz(vertex, momentum, bs.position());
778 
779  if(theSimPVPosition) {
780  dxyPVSim = TrackingParticleIP::dxy(vertex, momentum, *theSimPVPosition);
781  dzPVSim = TrackingParticleIP::dz(vertex, momentum, *theSimPVPosition);
782  }
783  }
784  //If the TrackingParticle is comics, get the momentum and vertex at PCA
785  else
786  {
787  momentumTP = std::get<TrackingParticle::Vector>(momVert);
788  vertexTP = std::get<TrackingParticle::Point>(momVert);
789  dxySim = TrackingParticleIP::dxy(vertexTP, momentumTP, bs.position());
790  dzSim = TrackingParticleIP::dz(vertexTP, momentumTP, bs.position());
791 
792  // Do dxy and dz vs. PV make any sense for cosmics? I guess not
793  }
794  //---------- THE PART ABOVE HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
795 
796  // in the coming lines, histos are filled using as input
797  // - momentumTP
798  // - vertexTP
799  // - dxySim
800  // - dzSim
801  if(!doSimTrackPlots_)
802  continue;
803 
804  // ##############################################
805  // fill RecoAssociated SimTracks' histograms
806  // ##############################################
807  const reco::Track *matchedTrackPointer = nullptr;
808  const reco::Track *matchedSecondTrackPointer = nullptr;
809  unsigned int selectsLoose = mvaCollections.size();
810  unsigned int selectsHP = mvaCollections.size();
811  if(simRecColl.find(tpr) != simRecColl.end()){
812  auto const & rt = simRecColl[tpr];
813  if (rt.size()!=0) {
814  ats++; //This counter counts the number of simTracks that have a recoTrack associated
815  // isRecoMatched = true; // UNUSED
816  matchedTrackPointer = rt.begin()->first.get();
817  if(rt.size() >= 2) {
818  matchedSecondTrackPointer = (rt.begin()+1)->first.get();
819  }
820  LogTrace("TrackValidator") << "TrackingParticle #" << st
821  << " with pt=" << sqrt(momentumTP.perp2())
822  << " associated with quality:" << rt.begin()->second <<"\n";
823 
824  if(doMVAPlots_) {
825  // for each MVA we need to take the value of the track
826  // with largest MVA value (for the cumulative histograms)
827  //
828  // also identify the first MVA that possibly selects any
829  // track matched to this TrackingParticle, separately
830  // for loose and highPurity qualities
831  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
832  const auto& mva = *(mvaCollections[imva]);
833  const auto& qual = *(qualityMaskCollections[imva]);
834 
835  auto iMatch = rt.begin();
836  float maxMva = mva[iMatch->first.key()];
837  for(; iMatch!=rt.end(); ++iMatch) {
838  auto itrk = iMatch->first.key();
839  maxMva = std::max(maxMva, mva[itrk]);
840 
841  if(selectsLoose >= imva && trackSelected(qual[itrk], reco::TrackBase::loose))
842  selectsLoose = imva;
843  if(selectsHP >= imva && trackSelected(qual[itrk], reco::TrackBase::highPurity))
844  selectsHP = imva;
845  }
846  mvaValues.push_back(maxMva);
847  }
848  }
849  }
850  }else{
851  LogTrace("TrackValidator")
852  << "TrackingParticle #" << st
853  << " with pt,eta,phi: "
854  << sqrt(momentumTP.perp2()) << " , "
855  << momentumTP.eta() << " , "
856  << momentumTP.phi() << " , "
857  << " NOT associated to any reco::Track" << "\n";
858  }
859 
860 
861 
862 
863  int nSimHits = tp.numberOfTrackerHits();
864  int nSimLayers = nLayers_tPCeff[tpr];
865  int nSimPixelLayers = nPixelLayers_tPCeff[tpr];
866  int nSimStripMonoAndStereoLayers = nStripMonoAndStereoLayers_tPCeff[tpr];
867  histoProducerAlgo_->fill_recoAssociated_simTrack_histos(w,tp,momentumTP,vertexTP,dxySim,dzSim,dxyPVSim,dzPVSim,nSimHits,nSimLayers,nSimPixelLayers,nSimStripMonoAndStereoLayers,matchedTrackPointer,puinfo.getPU_NumInteractions(), dR, thePVposition, theSimPVPosition, bs.position(), mvaValues, selectsLoose, selectsHP);
868  mvaValues.clear();
869 
870  if(matchedTrackPointer && matchedSecondTrackPointer) {
871  histoProducerAlgo_->fill_duplicate_histos(w, *matchedTrackPointer, *matchedSecondTrackPointer);
872  }
873 
874  if(doSummaryPlots_) {
875  if(dRtpSelector(tp)) {
876  h_simul_coll[ww]->Fill(www);
877  if (matchedTrackPointer) {
878  h_assoc_coll[ww]->Fill(www);
879  }
880  }
881  }
882 
883 
884 
885 
886  } // End for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
887 
888  // ##############################################
889  // fill recoTracks histograms (LOOP OVER TRACKS)
890  // ##############################################
891  if(!doRecoTrackPlots_)
892  continue;
893  LogTrace("TrackValidator") << "\n# of reco::Tracks with "
894  << label[www].process()<<":"
895  << label[www].label()<<":"
896  << label[www].instance()
897  << ": " << trackCollection.size() << "\n";
898 
899  int sat(0); //This counter counts the number of recoTracks that are associated to SimTracks from Signal only
900  int at(0); //This counter counts the number of recoTracks that are associated to SimTracks
901  int rT(0); //This counter counts the number of recoTracks in general
902  int seed_fit_failed = 0;
903  size_t n_selTrack_dr = 0;
904 
905  //calculate dR for tracks
906  const edm::View<Track> *trackCollectionDr = &trackCollection;
908  trackCollectionDr = trackCollectionForDrCalculation.product();
909  }
910  declareDynArray(float, trackCollection.size(), dR_trk);
911  trackDR(trackCollection, *trackCollectionDr, dR_trk);
912 
913  for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
914  auto track = trackCollection.refAt(i);
915  rT++;
916  if(trackFromSeedFitFailed(*track)) ++seed_fit_failed;
917  if((*dRTrackSelector)(*track)) ++n_selTrack_dr;
918 
919  bool isSigSimMatched(false);
920  bool isSimMatched(false);
921  bool isChargeMatched(true);
922  int numAssocRecoTracks = 0;
923  int nSimHits = 0;
924  double sharedFraction = 0.;
925 
926  auto tpFound = recSimColl.find(track);
927  isSimMatched = tpFound != recSimColl.end();
928  if (isSimMatched) {
929  const auto& tp = tpFound->val;
930  nSimHits = tp[0].first->numberOfTrackerHits();
931  sharedFraction = tp[0].second;
932  if (tp[0].first->charge() != track->charge()) isChargeMatched = false;
933  if(simRecColl.find(tp[0].first) != simRecColl.end()) numAssocRecoTracks = simRecColl[tp[0].first].size();
934  at++;
935  for (unsigned int tp_ite=0;tp_ite<tp.size();++tp_ite){
936  TrackingParticle trackpart = *(tp[tp_ite].first);
937  if ((trackpart.eventId().event() == 0) && (trackpart.eventId().bunchCrossing() == 0)){
938  isSigSimMatched = true;
939  sat++;
940  break;
941  }
942  }
943  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
944  << " associated with quality:" << tp.begin()->second <<"\n";
945  } else {
946  LogTrace("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
947  << " NOT associated to any TrackingParticle" << "\n";
948  }
949 
950  // set MVA values for this track
951  // take also the indices of first MVAs to select by loose and
952  // HP quality
953  unsigned int selectsLoose = mvaCollections.size();
954  unsigned int selectsHP = mvaCollections.size();
955  if(doMVAPlots_) {
956  for(size_t imva=0; imva<mvaCollections.size(); ++imva) {
957  const auto& mva = *(mvaCollections[imva]);
958  const auto& qual = *(qualityMaskCollections[imva]);
959  mvaValues.push_back(mva[i]);
960 
961  if(selectsLoose >= imva && trackSelected(qual[i], reco::TrackBase::loose))
962  selectsLoose = imva;
963  if(selectsHP >= imva && trackSelected(qual[i], reco::TrackBase::highPurity))
964  selectsHP = imva;
965  }
966  }
967 
968  double dR=dR_trk[i];
969  histoProducerAlgo_->fill_generic_recoTrack_histos(w,*track, ttopo, bs.position(), thePVposition, theSimPVPosition, isSimMatched,isSigSimMatched, isChargeMatched, numAssocRecoTracks, puinfo.getPU_NumInteractions(), nSimHits, sharedFraction, dR, mvaValues, selectsLoose, selectsHP);
970  mvaValues.clear();
971 
972  if(doSummaryPlots_) {
973  h_reco_coll[ww]->Fill(www);
974  if(isSimMatched) {
975  h_assoc2_coll[ww]->Fill(www);
976  if(numAssocRecoTracks>1) {
977  h_looper_coll[ww]->Fill(www);
978  }
979  if(!isSigSimMatched) {
980  h_pileup_coll[ww]->Fill(www);
981  }
982  }
983  }
984 
985  // dE/dx
986  if (dodEdxPlots_) histoProducerAlgo_->fill_dedx_recoTrack_histos(w,track, v_dEdx);
987 
988 
989  //Fill other histos
990  if (!isSimMatched) continue;
991 
992  histoProducerAlgo_->fill_simAssociated_recoTrack_histos(w,*track);
993 
994  /* TO BE FIXED LATER
995  if (associators[ww]=="trackAssociatorByChi2"){
996  //association chi2
997  double assocChi2 = -tp.begin()->second;//in association map is stored -chi2
998  h_assochi2[www]->Fill(assocChi2);
999  h_assochi2_prob[www]->Fill(TMath::Prob((assocChi2)*5,5));
1000  }
1001  else if (associators[ww]=="quickTrackAssociatorByHits"){
1002  double fraction = tp.begin()->second;
1003  h_assocFraction[www]->Fill(fraction);
1004  h_assocSharedHit[www]->Fill(fraction*track->numberOfValidHits());
1005  }
1006  */
1007 
1008 
1009  if(doResolutionPlots_[www]) {
1010  //Get tracking particle parameters at point of closest approach to the beamline
1011  TrackingParticleRef tpr = tpFound->val.begin()->first;
1012  TrackingParticle::Vector momentumTP = parametersDefinerTP->momentum(event,setup,tpr);
1013  TrackingParticle::Point vertexTP = parametersDefinerTP->vertex(event,setup,tpr);
1014  int chargeTP = tpr->charge();
1015 
1016  histoProducerAlgo_->fill_ResoAndPull_recoTrack_histos(w,momentumTP,vertexTP,chargeTP,
1017  *track,bs.position());
1018  }
1019 
1020 
1021  //TO BE FIXED
1022  //std::vector<PSimHit> simhits=tpr.get()->trackPSimHit(DetId::Tracker);
1023  //nrecHit_vs_nsimHit_rec2sim[w]->Fill(track->numberOfValidHits(), (int)(simhits.end()-simhits.begin() ));
1024 
1025  } // End of for(View<Track>::size_type i=0; i<trackCollection.size(); ++i){
1026  mvaCollections.clear();
1027  qualityMaskCollections.clear();
1028 
1029  histoProducerAlgo_->fill_trackBased_histos(w,at,rT, n_selTrack_dr, n_selTP_dr);
1030  // Fill seed-specific histograms
1031  if(doSeedPlots_) {
1032  histoProducerAlgo_->fill_seed_histos(www, seed_fit_failed, trackCollection.size());
1033  }
1034 
1035 
1036  LogTrace("TrackValidator") << "Collection " << www << "\n"
1037  << "Total Simulated (selected): " << n_selTP_dr << "\n"
1038  << "Total Reconstructed (selected): " << n_selTrack_dr << "\n"
1039  << "Total Reconstructed: " << rT << "\n"
1040  << "Total Associated (recoToSim): " << at << "\n"
1041  << "Total Fakes: " << rT-at << "\n";
1042  } // End of for (unsigned int www=0;www<label.size();www++){
1043  } //END of for (unsigned int ww=0;ww<associators.size();ww++){
1044 
1045 }
#define LogDebug(id)
size
Write out results.
void trackDR(const edm::View< reco::Track > &trackCollection, const edm::View< reco::Track > &trackCollectionDr, DynArray< float > &dR_trk) const
unsigned int size_type
Definition: View.h:90
edm::EDGetTokenT< SimHitTPAssociationProducer::SimHitTPAssociationList > _simHitTpMapTag
int event() const
get the contents of the subdetector field (should be protected?)
std::vector< MonitorElement * > h_reco_coll
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx2Tag
bool trackFromSeedFitFailed(const reco::Track &track)
std::vector< edm::EDGetTokenT< reco::SimToRecoCollection > > associatormapStRs
const double w
Definition: UKUtility.cc:23
CosmicTrackingParticleSelector cosmictpSelector
Vector momentum() const
spatial momentum vector
const_iterator end() const
last iterator over the map (read only)
std::vector< MonitorElement * > h_simul_coll
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
const_iterator find(const key_type &k) const
find element with specified reference key
std::vector< std::vector< std::tuple< edm::EDGetTokenT< MVACollection >, edm::EDGetTokenT< QualityMaskCollection > > > > mvaQualityCollectionTokens_
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNStripStereoLayersToken_
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNLayersToken_
size_type size() const
T_AssociationMap associationMapFilterValues(const T_AssociationMap &map, const T_RefVector &valueRefs)
TrackingParticleSelector dRtpSelector
std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associatorTokens
const bool doPlotsOnlyForTruePV_
edm::EDGetTokenT< reco::BeamSpot > bsSrc
const reco::Vertex::Point * getRecoPVPosition(const edm::Event &event, const edm::Handle< TrackingVertexCollection > &htv) const
void tpParametersAndSelection(const TrackingParticleRefVector &tPCeff, const ParametersDefinerForTP &parametersDefinerTP, const edm::Event &event, const edm::EventSetup &setup, const reco::BeamSpot &bs, std::vector< std::tuple< TrackingParticle::Vector, TrackingParticle::Point > > &momVert_tPCeff, std::vector< size_t > &selected_tPCeff) const
RefToBase< value_type > refAt(size_type i) const
math::XYZPointD Point
point in the space
std::vector< MonitorElement * > h_looper_coll
edm::EDGetTokenT< edm::ValueMap< unsigned int > > tpNPixelLayersToken_
math::XYZTLorentzVectorD LorentzVector
std::vector< edm::InputTag > associators
std::vector< MonitorElement * > h_pileup_coll
std::vector< MonitorElement * > h_assoc_coll
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > labelToken
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx1Tag
std::vector< edm::EDGetTokenT< reco::RecoToSimCollection > > associatormapRtSs
auto dz(const T_Vertex &vertex, const T_Momentum &momentum, const T_Point &point)
T sqrt(T t)
Definition: SSEVec.h:18
size_t tpDR(const TrackingParticleRefVector &tPCeff, const std::vector< size_t > &selected_tPCeff, DynArray< float > &dR_tPCeff) const
int bunchCrossing() const
get the detector field from this detid
std::vector< bool > doResolutionPlots_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< edm::View< reco::Track > > labelTokenForDrCalculation
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > label_pileupinfo
std::vector< edm::InputTag > label
void initEvent(edm::Handle< SimHitTPAssociationProducer::SimHitTPAssociationList > simHitsTPAssocToSet) const
virtual std::unique_ptr< ParametersDefinerForTP > clone() const
#define LogTrace(id)
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_fake_refvector
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< TrackingVertexCollection > label_tv
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
std::unique_ptr< RecoTrackSelectorBase > dRTrackSelector
const T & get() const
Definition: EventSetup.h:55
const bool ignoremissingtkcollection_
const int getPU_NumInteractions() const
const bool doPVAssociationPlots_
EncodedEventId eventId() const
Signal source, crossing number.
Point vertex() const
Parent vertex position.
const bool calculateDrSingleCollection_
fixed size matrix
auto dxy(const T_Vertex &vertex, const T_Momentum &momentum, const T_Point &point)
void push_back(const RefToBase< T > &)
std::vector< MonitorElement * > h_assoc2_coll
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
edm::EDGetTokenT< TrackingParticleCollection > label_tp_fake
Monte Carlo truth information used for tracking validation.
bool isUninitialized() const
Definition: EDGetToken.h:73
const Point & position() const
position
Definition: BeamSpot.h:62
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
math::XYZVectorD Vector
point in the space
edm::EDGetTokenT< TrackingParticleRefVector > label_tp_effic_refvector
edm::EDGetTokenT< TrackingParticleCollection > label_tp_effic
edm::Ref< TrackingParticleCollection > TrackingParticleRef
def move(src, dest)
Definition: eostools.py:510
const TrackingVertex::LorentzVector * getSimPVPosition(const edm::Handle< TrackingVertexCollection > &htv) const
const bool parametersDefinerIsCosmic_
void MultiTrackValidator::bookHistograms ( DQMStore::IBooker ibook,
edm::Run const &  ,
edm::EventSetup const &  setup 
)
overridevirtual

Method called to book the DQM histograms.

Implements DQMEDAnalyzer.

Definition at line 213 of file MultiTrackValidator.cc.

References patPFMETCorrections_cff::algo, trackingPlots::assoc, associators, DQMStore::IBooker::book1D(), DQMStore::IBooker::cd(), TrackerOfflineValidation_Dqm_cff::dirName, dirName_, dodEdxPlots_, doMVAPlots_, doPVAssociationPlots_, doRecoTrackPlots_, doResolutionPlots_, doSeedPlots_, doSimPlots_, doSimTrackPlots_, doSummaryPlots_, edm::RefVector< C, T, F >::empty(), Exception, h_assoc2_coll, h_assoc_coll, h_looper_coll, h_pileup_coll, h_reco_coll, h_simul_coll, edm::IndexSet::has(), histoProducerAlgo_, mps_fire::i, edm::RefVector< C, T, F >::id(), edm::IndexSet::insert(), edm::InputTag::instance(), edm::InputTag::label(), label, mvaQualityCollectionTokens_, edm::InputTag::process(), python.rootplot.root2matplotlib::replace(), edm::IndexSet::reserve(), DQMStore::IBooker::setCurrentFolder(), edm::RefVector< C, T, F >::size(), and findQualityFiles::size.

213  {
214 
215  const auto minColl = -0.5;
216  const auto maxColl = label.size()-0.5;
217  const auto nintColl = label.size();
218 
219  auto binLabels = [&](MonitorElement *me) {
220  TH1 *h = me->getTH1();
221  for(size_t i=0; i<label.size(); ++i) {
222  h->GetXaxis()->SetBinLabel(i+1, label[i].label().c_str());
223  }
224  return me;
225  };
226 
227  //Booking histograms concerning with simulated tracks
228  if(doSimPlots_) {
229  ibook.cd();
230  ibook.setCurrentFolder(dirName_ + "simulation");
231 
232  histoProducerAlgo_->bookSimHistos(ibook);
233 
234  ibook.cd();
235  ibook.setCurrentFolder(dirName_);
236  }
237 
238  for (unsigned int ww=0;ww<associators.size();ww++){
239  ibook.cd();
240  // FIXME: these need to be moved to a subdirectory whose name depends on the associator
241  ibook.setCurrentFolder(dirName_);
242 
243  if(doSummaryPlots_) {
244  if(doSimTrackPlots_) {
245  h_assoc_coll.push_back(binLabels( ibook.book1D("num_assoc(simToReco)_coll", "N of associated (simToReco) tracks vs track collection", nintColl, minColl, maxColl) ));
246  h_simul_coll.push_back(binLabels( ibook.book1D("num_simul_coll", "N of simulated tracks vs track collection", nintColl, minColl, maxColl) ));
247  }
248  if(doRecoTrackPlots_) {
249  h_reco_coll.push_back(binLabels( ibook.book1D("num_reco_coll", "N of reco track vs track collection", nintColl, minColl, maxColl) ));
250  h_assoc2_coll.push_back(binLabels( ibook.book1D("num_assoc(recoToSim)_coll", "N of associated (recoToSim) tracks vs track collection", nintColl, minColl, maxColl) ));
251  h_looper_coll.push_back(binLabels( ibook.book1D("num_duplicate_coll", "N of associated (recoToSim) looper tracks vs track collection", nintColl, minColl, maxColl) ));
252  h_pileup_coll.push_back(binLabels( ibook.book1D("num_pileup_coll", "N of associated (recoToSim) pileup tracks vs track collection", nintColl, minColl, maxColl) ));
253  }
254  }
255 
256  for (unsigned int www=0;www<label.size();www++){
257  ibook.cd();
258  InputTag algo = label[www];
259  string dirName=dirName_;
260  if (algo.process()!="")
261  dirName+=algo.process()+"_";
262  if(algo.label()!="")
263  dirName+=algo.label()+"_";
264  if(algo.instance()!="")
265  dirName+=algo.instance()+"_";
266  if (dirName.find("Tracks")<dirName.length()){
267  dirName.replace(dirName.find("Tracks"),6,"");
268  }
269  string assoc= associators[ww].label();
270  if (assoc.find("Track")<assoc.length()){
271  assoc.replace(assoc.find("Track"),5,"");
272  }
273  dirName+=assoc;
274  std::replace(dirName.begin(), dirName.end(), ':', '_');
275 
276  ibook.setCurrentFolder(dirName.c_str());
277 
278  const bool doResolutionPlots = doResolutionPlots_[www];
279 
280  if(doSimTrackPlots_) {
281  histoProducerAlgo_->bookSimTrackHistos(ibook, doResolutionPlots);
282  if(doPVAssociationPlots_) histoProducerAlgo_->bookSimTrackPVAssociationHistos(ibook);
283  }
284 
285  //Booking histograms concerning with reconstructed tracks
286  if(doRecoTrackPlots_) {
287  histoProducerAlgo_->bookRecoHistos(ibook, doResolutionPlots);
288  if (dodEdxPlots_) histoProducerAlgo_->bookRecodEdxHistos(ibook);
289  if (doPVAssociationPlots_) histoProducerAlgo_->bookRecoPVAssociationHistos(ibook);
290  if (doMVAPlots_) histoProducerAlgo_->bookMVAHistos(ibook, mvaQualityCollectionTokens_[www].size());
291  }
292 
293  if(doSeedPlots_) {
294  histoProducerAlgo_->bookSeedHistos(ibook);
295  }
296  }//end loop www
297  }// end loop ww
298 }
size
Write out results.
std::vector< MonitorElement * > h_reco_coll
std::vector< MonitorElement * > h_simul_coll
void cd(void)
Definition: DQMStore.cc:269
std::vector< std::vector< std::tuple< edm::EDGetTokenT< MVACollection >, edm::EDGetTokenT< QualityMaskCollection > > > > mvaQualityCollectionTokens_
def replace(string, replacements)
std::vector< MonitorElement * > h_looper_coll
std::vector< edm::InputTag > associators
std::vector< MonitorElement * > h_pileup_coll
std::vector< MonitorElement * > h_assoc_coll
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
std::vector< bool > doResolutionPlots_
std::vector< edm::InputTag > label
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const bool doPVAssociationPlots_
std::string const & label() const
Definition: InputTag.h:36
std::string const & process() const
Definition: InputTag.h:40
std::vector< MonitorElement * > h_assoc2_coll
std::string const & instance() const
Definition: InputTag.h:37
const reco::Vertex::Point * MultiTrackValidator::getRecoPVPosition ( const edm::Event event,
const edm::Handle< TrackingVertexCollection > &  htv 
) const
private

Definition at line 338 of file MultiTrackValidator.cc.

References reco::VertexToTrackingVertexAssociator::associateRecoToSim(), EncodedEventId::bunchCrossing(), EncodedEventId::event(), TrackingVertex::eventId(), recoVertexToken_, and vertexAssociatorToken_.

Referenced by analyze().

338  {
340  event.getByToken(recoVertexToken_, hvertex);
341 
343  event.getByToken(vertexAssociatorToken_, hvassociator);
344 
345  auto v_r2s = hvassociator->associateRecoToSim(hvertex, htv);
346  auto pvPtr = hvertex->refAt(0);
347  if(pvPtr->isFake() || pvPtr->ndof() < 0) // skip junk vertices
348  return nullptr;
349 
350  auto pvFound = v_r2s.find(pvPtr);
351  if(pvFound == v_r2s.end())
352  return nullptr;
353 
354  for(const auto& vertexRefQuality: pvFound->val) {
355  const TrackingVertex& tv = *(vertexRefQuality.first);
356  if(tv.eventId().event() == 0 && tv.eventId().bunchCrossing() == 0) {
357  return &(pvPtr->position());
358  }
359  }
360 
361  return nullptr;
362 }
int event() const
get the contents of the subdetector field (should be protected?)
int bunchCrossing() const
get the detector field from this detid
edm::EDGetTokenT< reco::VertexToTrackingVertexAssociator > vertexAssociatorToken_
edm::EDGetTokenT< edm::View< reco::Vertex > > recoVertexToken_
const EncodedEventId & eventId() const
reco::VertexRecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Vertex > > &vCH, const edm::Handle< TrackingVertexCollection > &tVCH) const
compare reco to sim the handle of reco::Vertex and TrackingVertex collections
const TrackingVertex::LorentzVector * MultiTrackValidator::getSimPVPosition ( const edm::Handle< TrackingVertexCollection > &  htv) const
private

Definition at line 329 of file MultiTrackValidator.cc.

Referenced by analyze().

329  {
330  for(const auto& simV: *htv) {
331  if(simV.eventId().bunchCrossing() != 0) continue; // remove OOTPU
332  if(simV.eventId().event() != 0) continue; // pick the PV of hard scatter
333  return &(simV.position());
334  }
335  return nullptr;
336 }
size_t MultiTrackValidator::tpDR ( const TrackingParticleRefVector tPCeff,
const std::vector< size_t > &  selected_tPCeff,
DynArray< float > &  dR_tPCeff 
) const
private

Definition at line 421 of file MultiTrackValidator.cc.

References reco::deltaR2(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, dRtpSelector, PVValHelper::eta, SiStripPI::max, AlCaHLTBitMon_ParallelJobs::p, phi, edm::RefVector< C, T, F >::size(), and mathSSE::sqrt().

Referenced by analyze().

423  {
424  float etaL[tPCeff.size()], phiL[tPCeff.size()];
425  size_t n_selTP_dr = 0;
426  for(size_t iTP: selected_tPCeff) {
427  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
428  auto const& tp2 = *(tPCeff[iTP]);
429  auto && p = tp2.momentum();
430  etaL[iTP] = etaFromXYZ(p.x(),p.y(),p.z());
431  phiL[iTP] = atan2f(p.y(),p.x());
432  }
433  for(size_t iTP1: selected_tPCeff) {
434  auto const& tp = *(tPCeff[iTP1]);
436  if(dRtpSelector(tp)) {//only for those needed for efficiency!
437  ++n_selTP_dr;
438  float eta = etaL[iTP1];
439  float phi = phiL[iTP1];
440  for(size_t iTP2: selected_tPCeff) {
441  //calculare dR wrt inclusive collection (also with PU, low pT, displaced)
442  if (iTP1==iTP2) {continue;}
443  auto dR_tmp = reco::deltaR2(eta, phi, etaL[iTP2], phiL[iTP2]);
444  if (dR_tmp<dR) dR=dR_tmp;
445  } // ttp2 (iTP)
446  }
447  dR_tPCeff[iTP1] = std::sqrt(dR);
448  } // tp
449  return n_selTP_dr;
450 }
TrackingParticleSelector dRtpSelector
T sqrt(T t)
Definition: SSEVec.h:18
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
void MultiTrackValidator::tpParametersAndSelection ( const TrackingParticleRefVector tPCeff,
const ParametersDefinerForTP parametersDefinerTP,
const edm::Event event,
const edm::EventSetup setup,
const reco::BeamSpot bs,
std::vector< std::tuple< TrackingParticle::Vector, TrackingParticle::Point > > &  momVert_tPCeff,
std::vector< size_t > &  selected_tPCeff 
) const
private

Definition at line 364 of file MultiTrackValidator.cc.

References EncodedEventId::bunchCrossing(), cosmictpSelector, doSimPlots_, TrackingParticle::eventId(), histoProducerAlgo_, ParametersDefinerForTP::momentum(), TrackingParticle::momentum(), parametersDefinerIsCosmic_, edm::RefVector< C, T, F >::size(), tpSelector, ParametersDefinerForTP::vertex(), and TrackingParticle::vertex().

Referenced by analyze().

369  {
370  selected_tPCeff.reserve(tPCeff.size());
371  momVert_tPCeff.reserve(tPCeff.size());
372  int nIntimeTPs = 0;
374  for(size_t j=0; j<tPCeff.size(); ++j) {
375  const TrackingParticleRef& tpr = tPCeff[j];
376 
377  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
378  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
379  if(doSimPlots_) {
380  histoProducerAlgo_->fill_generic_simTrack_histos(momentum, vertex, tpr->eventId().bunchCrossing());
381  }
382  if(tpr->eventId().bunchCrossing() == 0)
383  ++nIntimeTPs;
384 
385  if(cosmictpSelector(tpr,&bs,event,setup)) {
386  selected_tPCeff.push_back(j);
387  momVert_tPCeff.emplace_back(momentum, vertex);
388  }
389  }
390  }
391  else {
392  size_t j=0;
393  for(auto const& tpr: tPCeff) {
394  const TrackingParticle& tp = *tpr;
395 
396  // TODO: do we want to fill these from all TPs that include IT
397  // and OOT (as below), or limit to IT+OOT TPs passing tpSelector
398  // (as it was before)? The latter would require another instance
399  // of tpSelector with intimeOnly=False.
400  if(doSimPlots_) {
401  histoProducerAlgo_->fill_generic_simTrack_histos(tp.momentum(), tp.vertex(), tp.eventId().bunchCrossing());
402  }
403  if(tp.eventId().bunchCrossing() == 0)
404  ++nIntimeTPs;
405 
406  if(tpSelector(tp)) {
407  selected_tPCeff.push_back(j);
408  TrackingParticle::Vector momentum = parametersDefinerTP.momentum(event,setup,tpr);
409  TrackingParticle::Point vertex = parametersDefinerTP.vertex(event,setup,tpr);
410  momVert_tPCeff.emplace_back(momentum, vertex);
411  }
412  ++j;
413  }
414  }
415  if(doSimPlots_) {
416  histoProducerAlgo_->fill_simTrackBased_histos(nIntimeTPs);
417  }
418 }
CosmicTrackingParticleSelector cosmictpSelector
Vector momentum() const
spatial momentum vector
math::XYZPointD Point
point in the space
int bunchCrossing() const
get the detector field from this detid
TrackingParticleSelector tpSelector
virtual TrackingParticle::Vector momentum(const edm::Event &iEvent, const edm::EventSetup &iSetup, const Charge ch, const Point &vtx, const LorentzVector &lv) const
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
EncodedEventId eventId() const
Signal source, crossing number.
Point vertex() const
Parent vertex position.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
Monte Carlo truth information used for tracking validation.
math::XYZVectorD Vector
point in the space
virtual TrackingParticle::Point vertex(const edm::Event &iEvent, const edm::EventSetup &iSetup, const Charge ch, const Point &vtx, const LorentzVector &lv) const
const bool parametersDefinerIsCosmic_
void MultiTrackValidator::trackDR ( const edm::View< reco::Track > &  trackCollection,
const edm::View< reco::Track > &  trackCollectionDr,
DynArray< float > &  dR_trk 
) const
private

Definition at line 452 of file MultiTrackValidator.cc.

References reco::deltaR2(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, PVValHelper::eta, mps_fire::i, SiStripPI::max, min(), AlCaHLTBitMon_ParallelJobs::p, phi, edm::View< T >::size(), mathSSE::sqrt(), HiIsolationCommonParameters_cff::track, and trackFromSeedFitFailed().

Referenced by analyze().

452  {
453  int i=0;
454  float etaL[trackCollectionDr.size()];
455  float phiL[trackCollectionDr.size()];
456  bool validL[trackCollectionDr.size()];
457  for (auto const & track2 : trackCollectionDr) {
458  auto && p = track2.momentum();
459  etaL[i] = etaFromXYZ(p.x(),p.y(),p.z());
460  phiL[i] = atan2f(p.y(),p.x());
461  validL[i] = !trackFromSeedFitFailed(track2);
462  ++i;
463  }
464  for(View<reco::Track>::size_type i=0; i<trackCollection.size(); ++i){
465  auto const & track = trackCollection[i];
468  auto && p = track.momentum();
469  float eta = etaFromXYZ(p.x(),p.y(),p.z());
470  float phi = atan2f(p.y(),p.x());
471  for(View<reco::Track>::size_type j=0; j<trackCollectionDr.size(); ++j){
472  if(!validL[j]) continue;
473  auto dR_tmp = reco::deltaR2(eta, phi, etaL[j], phiL[j]);
474  if ( (dR_tmp<dR) & (dR_tmp>std::numeric_limits<float>::min())) dR=dR_tmp;
475  }
476  }
477  dR_trk[i] = std::sqrt(dR);
478  }
479 }
unsigned int size_type
Definition: View.h:90
bool trackFromSeedFitFailed(const reco::Track &track)
size_type size() const
T sqrt(T t)
Definition: SSEVec.h:18
T min(T a, T b)
Definition: MathUtil.h:58
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36

Member Data Documentation

edm::EDGetTokenT<SimHitTPAssociationProducer::SimHitTPAssociationList> MultiTrackValidator::_simHitTpMapTag
private

Definition at line 128 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<edm::EDGetTokenT<reco::RecoToSimCollection> > MultiTrackValidator::associatormapRtSs
private

Definition at line 106 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<edm::EDGetTokenT<reco::SimToRecoCollection> > MultiTrackValidator::associatormapStRs
private

Definition at line 105 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<edm::InputTag> MultiTrackValidator::associators
protected
std::vector<edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator> > MultiTrackValidator::associatorTokens
private

Definition at line 104 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

edm::EDGetTokenT<reco::BeamSpot> MultiTrackValidator::bsSrc
protected
const bool MultiTrackValidator::calculateDrSingleCollection_
protected

Definition at line 73 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

CosmicTrackingParticleSelector MultiTrackValidator::cosmictpSelector
private

Definition at line 124 of file MultiTrackValidator.h.

Referenced by analyze(), MultiTrackValidator(), and tpParametersAndSelection().

std::string MultiTrackValidator::dirName_
private

Definition at line 117 of file MultiTrackValidator.h.

Referenced by bookHistograms(), and MultiTrackValidator().

const bool MultiTrackValidator::dodEdxPlots_
protected
const bool MultiTrackValidator::doMVAPlots_
protected

Definition at line 82 of file MultiTrackValidator.h.

Referenced by analyze(), bookHistograms(), and MultiTrackValidator().

const bool MultiTrackValidator::doPlotsOnlyForTruePV_
protected

Definition at line 74 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

const bool MultiTrackValidator::doPVAssociationPlots_
protected

Definition at line 80 of file MultiTrackValidator.h.

Referenced by analyze(), bookHistograms(), and MultiTrackValidator().

const bool MultiTrackValidator::doRecoTrackPlots_
protected
std::vector<bool> MultiTrackValidator::doResolutionPlots_
protected

Definition at line 84 of file MultiTrackValidator.h.

Referenced by analyze(), bookHistograms(), and MultiTrackValidator().

const bool MultiTrackValidator::doSeedPlots_
protected

Definition at line 81 of file MultiTrackValidator.h.

Referenced by analyze(), bookHistograms(), and MultiTrackValidator().

const bool MultiTrackValidator::doSimPlots_
protected
const bool MultiTrackValidator::doSimTrackPlots_
protected
const bool MultiTrackValidator::doSummaryPlots_
protected

Definition at line 75 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

TrackingParticleSelector MultiTrackValidator::dRtpSelector
private

Definition at line 125 of file MultiTrackValidator.h.

Referenced by analyze(), MultiTrackValidator(), and tpDR().

std::unique_ptr<RecoTrackSelectorBase> MultiTrackValidator::dRTrackSelector
private

Definition at line 126 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<MonitorElement *> MultiTrackValidator::h_assoc2_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> MultiTrackValidator::h_assoc_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> MultiTrackValidator::h_looper_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> MultiTrackValidator::h_pileup_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> MultiTrackValidator::h_reco_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> MultiTrackValidator::h_simul_coll
private

Definition at line 133 of file MultiTrackValidator.h.

Referenced by analyze(), and bookHistograms().

std::unique_ptr<MTVHistoProducerAlgoForTracker> MultiTrackValidator::histoProducerAlgo_
protected
const bool MultiTrackValidator::ignoremissingtkcollection_
protected

Definition at line 71 of file MultiTrackValidator.h.

Referenced by MultiTrackValidatorGenPs::analyze(), and analyze().

std::vector<edm::InputTag> MultiTrackValidator::label
protected
edm::EDGetTokenT<std::vector<PileupSummaryInfo> > MultiTrackValidator::label_pileupinfo
protected
edm::EDGetTokenT<TrackingParticleCollection> MultiTrackValidator::label_tp_effic
protected
edm::EDGetTokenT<TrackingParticleRefVector> MultiTrackValidator::label_tp_effic_refvector
protected

Definition at line 53 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

edm::EDGetTokenT<TrackingParticleCollection> MultiTrackValidator::label_tp_fake
protected
edm::EDGetTokenT<TrackingParticleRefVector> MultiTrackValidator::label_tp_fake_refvector
protected

Definition at line 54 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

edm::EDGetTokenT<TrackingVertexCollection> MultiTrackValidator::label_tv
protected

Definition at line 55 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<edm::EDGetTokenT<edm::View<reco::Track> > > MultiTrackValidator::labelToken
protected
edm::EDGetTokenT<edm::View<reco::Track> > MultiTrackValidator::labelTokenForDrCalculation
private

Definition at line 129 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

std::vector<edm::EDGetTokenT<edm::View<TrajectorySeed> > > MultiTrackValidator::labelTokenSeed
protected

Definition at line 62 of file MultiTrackValidator.h.

edm::EDGetTokenT<edm::ValueMap<reco::DeDxData> > MultiTrackValidator::m_dEdx1Tag
protected
edm::EDGetTokenT<edm::ValueMap<reco::DeDxData> > MultiTrackValidator::m_dEdx2Tag
protected
std::vector<std::vector<std::tuple<edm::EDGetTokenT<MVACollection>, edm::EDGetTokenT<QualityMaskCollection> > > > MultiTrackValidator::mvaQualityCollectionTokens_
private

Definition at line 115 of file MultiTrackValidator.h.

Referenced by analyze(), bookHistograms(), and MultiTrackValidator().

std::string MultiTrackValidator::parametersDefiner
protected

Definition at line 68 of file MultiTrackValidator.h.

Referenced by MultiTrackValidatorGenPs::analyze(), and analyze().

const bool MultiTrackValidator::parametersDefinerIsCosmic_
protected
edm::EDGetTokenT<edm::View<reco::Vertex> > MultiTrackValidator::recoVertexToken_
private

Definition at line 130 of file MultiTrackValidator.h.

Referenced by getRecoPVPosition(), and MultiTrackValidator().

std::vector<edm::EDGetTokenT<std::vector<PSimHit> > > MultiTrackValidator::simHitTokens_
protected

Definition at line 58 of file MultiTrackValidator.h.

Referenced by MultiTrackValidator().

const double MultiTrackValidator::simPVMaxZ_
private

Definition at line 120 of file MultiTrackValidator.h.

Referenced by analyze().

edm::EDGetTokenT<edm::ValueMap<unsigned int> > MultiTrackValidator::tpNLayersToken_
private

Definition at line 108 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

edm::EDGetTokenT<edm::ValueMap<unsigned int> > MultiTrackValidator::tpNPixelLayersToken_
private

Definition at line 109 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

edm::EDGetTokenT<edm::ValueMap<unsigned int> > MultiTrackValidator::tpNStripStereoLayersToken_
private

Definition at line 110 of file MultiTrackValidator.h.

Referenced by analyze(), and MultiTrackValidator().

TrackingParticleSelector MultiTrackValidator::tpSelector
private

Definition at line 123 of file MultiTrackValidator.h.

Referenced by MultiTrackValidator(), and tpParametersAndSelection().

const bool MultiTrackValidator::useAssociators_
protected
bool MultiTrackValidator::useGsf
private

Definition at line 119 of file MultiTrackValidator.h.

Referenced by MultiTrackValidator().

edm::EDGetTokenT<reco::VertexToTrackingVertexAssociator> MultiTrackValidator::vertexAssociatorToken_
private

Definition at line 131 of file MultiTrackValidator.h.

Referenced by getRecoPVPosition(), and MultiTrackValidator().