MuonAnalysis/MuonAssociators package documentation

Introduction

This twikipage will document the tools and cmssw modules in the MuonAnalysis/MuonAssociators, developed to associate muons or other muon-like objects to other information (trigger muons, generator level muons, ...).

The algorithms themselves can be used from any CMSSW module and carry no dependency on PAT or PhysicsTools (except the StringParser, which is in PhysicsTools in release 2.2.X)

The cmssw modules under MuonAnalysis/MuonAssociators/plugins are written with the physics analysis usecase in mind, and to be smoothly interfaced to PAT (even if they can be used without PAT)

The code is available for:

  • CMSSW_3_8_X and later: tag V01-12-01 (tested on 3_8_7, 3_9_7)
  • CMSSW_3_5_X and later: tag V01-07-05 (tested on 3_5_8, 3_6_1_patch2)
  • CMSSW_3_4_2 and later: tag V01-03-00 (tested on 3_4_2)
  • CMSSW_3_1_X and later: tag V01-02-00 (tested on 3_1_4 , 3_4_1; probably V01-03-00 works too)
  • CMSSW_2_2_X : V00-01-02 (tested on 2_2_10 )

Trigger-related code

Propagator to Muon Station 2

This object propagates the track or 4-momentum of a Candidate object up to the second muon station, so that it can be matched to the L1

The interface provided is TrajectoryStateOnSurface extrapolate(something) ; it will return a valid trajectory state there, or a null state if the propagation fails. The something can be a reco::Candidate, reco::Track or FreeTrajectoryState

Note that the thing to compare with the L1 eta and phi are the onse from the state globalPosition, not globalMomentum!

Example usage:

/// in the beginRun or at the beginning of the event
   propagatorToMuon.init(iSetup);
   /// then
   const reco::Muon &mu = ...
   TrajectoryStateOnSurface = propagatorToMuon.extrapolate(mu);
   if (prop.isValid()) {
      std::cout << "Propagation suceeded, eta = " << prop.globalPosition().eta() << ", phi = " << prop.globalPosition().phi() << std::endl;
   }

Conditions: At the beginning of each run (or each event), you must call the method void init(const edm::EventSetup & iSetup) so that the matcher can take from the EventSetup the geometry and the propagators. The following includes might be need in your cfg file, in addition to the standard ones like conditions, global tag, geometry and magnetic field:

process.load("TrackPropagation.SteppingHelixPropagator.SteppingHelixPropagatorAny_cfi")
process.load("TrackPropagation.SteppingHelixPropagator.SteppingHelixPropagatorAlong_cfi")
process.load("TrackPropagation.SteppingHelixPropagator.SteppingHelixPropagatorOpposite_cfi")
process.load("RecoMuon.DetLayers.muonDetLayerGeometry_cfi")

The propagator takes the following configuration parameters:

name type required notes
useTrack cms.string yes one of none, tracker, muon, global. It defines the type of the track that should be extracted from the Candidate object.
none means to just make a fake track with the candidate's momentum and charge.
Note that objects of type RecoChanrgedCandidate , like the ones used by L2 or L3 muon reco, will believe they have a tracker track, not a muon or global one.
useState cms.string sometimes Define which state to start from: atVertex (aka at PCA), innermost, outermost.
Required if useTrack is not none
Note that anything except atVertex will require the reco::TrackExtra to be available.
useSimpleGeometry cms.bool yes If true, the muon station geometry is taken directly from the DetLayers (i.e. one cylinder and two disks),
if false, a further propagation to the individual chambers will be attempted.
The Recommended value is true.
cosmicPropagationHypothesis cms.bool no If present and set to true, tracks will be considered cosmics, and the propagation direction will be computed accordingly.
Warning, important This mode has not been carefully tested.

Code: interface/PropagateToMuon.h , src/PropagateToMuon.cc

Note: before tag V00-01-00 this object was one with the L1MuonMatcherAlgo.

Associator Reco to L1 muons

This associator maps reconstructed objects to L1 muon candidates, by using the propagator to muon station 2 described before, and performing a match by deltaR and deltaPhi there.

Algorithm

The algorithm object takes care of the extrapolation and of the matching.

The three main methods provided to the user are, from lower level to higher level:

  • bool match(something, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated)
    Propagate the something to the muon station 2 and check if it matches with the object; the something can be a reco::Candidate or reco::Track . If the propagation succeeds, the propagated state will be updated as for the extrapolate method described above.
    If the matching succeeds, the deltaR and deltaPhi parameters will also be updated.
    If a preselection cut was specified in the configuration of the matcher, and the L1 fails the cut, the match will always be considered a failure (but not the propagation)
  • int match(something, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated)
    Same as above, but considers many L1 candidates and return the index of the best match or -1 if no match is found.

For backwards compatibility, it exposes the TrajectoryStateOnSurface extrapolate(something) from the propagator object.

Example usage:

edm::Handle<std::vector<l1extra::L1MuonParticle> > l1s;
   iEvent.getByLabel(l1s_, l1s);
   const reco::Muon &mu = ...
   float deltaR, deltaPhi; TrajectoryStateOnSurface prop;
   int match = matcher_.match(mu, *l1s, deltaR, deltaPhi, prop);
   if (prop.isValid()) {
      std::cout << "Propagation suceeded, eta = " << prop.globalPosition().eta() << ", phi = " << prop.globalPosition().phi() << std::endl;
   }
   if (match != -1) {
      std::cout << "Matched. L1 pt code is  " << (*l1s)[match].pt() << ", deltaPhi = " << deltaPhi << std::endl;
   }

Conditions: Just like for the PropagateToMuon algorithm object:

  • at the beginning of each run (or each event), you must call the method void init(const edm::EventSetup & iSetup) so that the matcher can take from the EventSetup the geometry and the propagators.
  • you need to have in your cfg file the necessary includes to define the needed EventSetup services (magnetic field, geometry, propagators)

The matcher takes all the parameters of the PropagateToMuon object, plus the following configuration parameters:

name type required notes
useTrack cms.string yes one of none, tracker, muon, global. It defines the type of the track that should be extracted from the candidate. none means to use the candidate's momentum vector and charge
preselection cms.string no A preselection cut to be applied to the L1MuonParticle objects. If not specified, no preselection will be done.
maxDeltaPhi cms.double noA Maximum allowed deltaPhi between the L1 candidate and the extrapolated reco object
maxDeltaR cms.double yes Maximum allowed deltaPhi between the L1 candidate and the extrapolated reco object
sortByDeltaPhi cms.bool no Sort matches by deltaPhi instead of deltaR (if not specified, it will be considered false)

(A) it was required in tags older than V01-03-00.

Code: interface/L1MuonMatcherAlgo.h , src/L1MuonMatcherAlgo.cc

Warning, important Known issues: Warning, important

  • Warning, important Before tag V00-00-04, the matching by deltaPhi is wrong (it's doing deltaPhi < cut instead of fabs(deltaPhi) < cut)

CMSSW Module

The module runs the matcher and produces some Associations and ValueMaps that can be used from an EDAnalyzer or copied into a pat::Muon, pat::GenericParticle.

Code: plugins/L1MuonMatcher.cc Config: muonL1Match_cfi.py

A full list of the meaningful configurables is

name type required notes
src cms.InputTag yes the reco objects, something that can be read through =View
matched cms.InputTag yes the collection of L1MuonParticle
preselection cms.string no A preselection cut to be applied to the L1MuonParticle objects. If not specified, no preselection will be done.
useTrack cms.string yes one of none, tracker, muon, global. It defines the type of the track that should be extracted from the Candidate object.
none means to just make a fake track with the candidate's momentum and charge.
Note that objects of type RecoChanrgedCandidate , like the ones used by L2 or L3 muon reco, will believe they have a tracker track, not a muon or global one.
useState cms.string sometimes Define which state to start from: atVertex (aka at PCA), innermost, outermost.
Required if useTrack is not none
Note that anything except atVertex will require the reco::TrackExtra to be available.
useSimpleGeometry cms.bool yes If true, the muon station geometry is taken directly from the DetLayers (i.e. one cylinder and two disks),
if false, a further propagation to the individual chambers will be attempted.
The Recommended value is true.
maxDeltaR cms.double yes Maximum allowed deltaR between the L1 candidate and the extrapolated reco object
maxDeltaPhi cms.double no If you want, you can specify also a maximum deltaPhi
sortByDeltaPhi cms.bool no Sort matches by deltaPhi instead of deltaR (if not specified, it will be considered false)
setPropLabel cms.string yes The fake filter label to assign to the pat::TriggerObjectStandAlone made to wrap the 4-momentum of the reco object propagated to the second muon station
setL1Label cms.string yes The fake filter label to assign to the pat::TriggerObjectStandAlone made to wrap the L1MuonParticle objects passing the preselection.

Basic example: run the L1 matcher with the default settings, and save the match result in the pat::Muon among the standard PAT trigger matches.

process.load("MuonAnalysis.MuonAssociators.muonL1Match_cfi")
process.patDefaultSequence.replace(process.muonTrigMatchHLT1MuonNonIso, process.muonTrigMatchHLT1MuonNonIso + process.muonL1Match)
process.allLayer1Muons.trigPrimMatch += [ cms.InputTag("muonL1Match") ]

Within your analyzer code, you can then check if the muon matches simply by

std::vector<pat::TriggerObjectStandAlone> match = muon.triggerMatchesByFilter("l1");
if (match.empty()) { 
   // sorry, no match 
} else { 
   float l1PtCode = match[0].pt(); 
}

Note that the pt, eta and phi are defined as for L1 candidates: eta and phi are positions at muon station 2, pt is at vertex but it's not really an estimate of the track pt but rather a lower bound on the muon pt with a given CL (i.e. 95% or 90% of the muons with pt at vertex > 5 GeV will have a L1 pt code > 5 GeV)

Warning, important Notes:

  • If the requested track for extrapolation is null, the muon will be considered unmatched. This default configuration in muonL1Match_cfi uses the tracker track, so you shouldn't expect matches for muon without a tracker track.
    You can also use the standalone muon track, the global track or the bare momentum and charge of the Candidate (the latter might be needed if you're trying to match an object that doesn't have a reconstructed track). You can also choose to start from the state at vertex, the innermost measurement or the outermost measurement on the track; anyway, take into account that for tracker tracks only the state at vertex is saved in the AOD data tier.
  • In the default configuration, l1extraParticles are produced only for the L1A bunch crossing. Anyway, you can re-make l1extraParticles for muons also in the previous and next bunch crossing, as in the example config file test/L1MuonMatcher/testAnyBX.py.
  • The special muon propagation feature must be enabled explicitly by adding the parameter cosmicPropagationHypothesis = cms.bool(True) as shown in test/L1MuonMatcher/testCRAFTCosmic.py. This feature has not been well tested yet, so while it looks like there are no obvious problems (e.g. it matches both top and bottom muons) still there might be hidden bugs or inefficiencies.

Advanced example: here we take also all the extra information about the L1 candidate, and put it into the CMS.UserData fields of the pat::Muon; we also keep the information about the propagated RECO track.

process.load("MuonAnalysis.MuonAssociators.muonL1Match_cfi")
process.muonL1Match.preselection = cms.string("")
process.muonL1Match.writeExtraInfo = True
process.patDefaultSequence.replace(process.muonTrigMatchHLT1MuonNonIso, process.muonL1Match)
process.allLayer1Muons.trigPrimMatch = cms.VInputTag(
    cms.InputTag("muonL1Match"),  # the match with the real L1 candidate
    cms.InputTag("muonL1Match","propagatedReco"), # the 4-vector from the propagation
)
process.allLayer1Muons.userData.userInts.src = cms.VInputTag(  # all these will be -999 in case of no match
    cms.InputTag("muonL1Match", "bx"),
    cms.InputTag("muonL1Match", "quality"),
    cms.InputTag("muonL1Match", "isolated"),
)
process.allLayer1Muons.userData.userFloats.src = cms.VInputTag( # this will be 999 in case of no match
    cms.InputTag("muonL1Match", "deltaR")
)
In this case, you can e.g. compare the phi of the propagated track and of the L1, and check the quality value of the trigger candidate
std::vector<pat::TriggerObjectStandAlone> prop = muon.triggerMatchesByFilter("propagated");
std::vector<pat::TriggerObjectStandAlone> match = muon.triggerMatchesByFilter("l1");
if (prop.empty()) { 
   // sorry, it didn't reach station 2
} else if (match.empty()) {
   // sorry, it didn't match
} else { 
   float phiDiff = deltaPhi( prop[0].phi(), match[0].phi() ); 
   int quality = muon.userInt("muonL1Match:quality");
}

In case you don't want to use a pat::Muon or pat::!GenericParticle to collapse this information, you can just run the matcher and then access directly the ValueMaps:

Handle<View<reco::Muon> > muons; 
iEvent.getByLabel(muonTag_, muons);

Handle<Association<pat::TriggerObjectStandAlone> > patmatches; // matches in PAT format
iEvent.getByLabel(patmatcher_, patmatches);

Handle<ValueMap<reco::CandidatePtr> > rawmatches; // matches to the l1extra directly
iEvent.getByLabel(rawmatcher_, rawmatches);

Handle<ValueMap<int> > qualities;
iEvent.getByLabel(edm::InputTag(matcher_.label(), "quality"), qualities);


for (size_t i = 0, n = muons->size(); i < n; ++i) {
     Ptr<reco::Muon> muPtr = muons->ptrAt(i);
     // Access data through the PAT primitive object and the ValueMaps
     pat::TriggerObjectStandAloneRef patmatch = (*patmatches)[muPtr];
     if (patmatch.isNonnull()) {
          float l1pt = patmatch->pt();
          int quality = (*qualities)[muPtr];
     } 
     // or, if you prefer, access the GMT Candidate directly
     reco::CandidatePtr rawmatchptr = (*rawmatches)[muPtr];
     if (rawmatchptr.isNonnull()) {
          const l1extra::L1MuonParticle  &l1p = dynamic_cast<const l1extra::L1MuonParticle &>(*rawmatchptr);
          const L1MuGMTExtendedCand &gmtCand = l1p.gmtMuonCand();
          bool rpc = gmtCand.isRPC();
          int bx= gmtCand.bx();
     } 
}

CMSSW Module for HLT paths which are L1 passthroughs (PAT-only)

When using PAT Trigger matching tools, the 4-vectors associated to HLT objects are converted into pat::TriggerObjectStandAlone objects and then are matched by δR to the reco objects. This doesn't work that well for HLT paths which are just passthroughs of L1 muons, because they will have the (η, φ) defined at the muon station 2 instead of the vertex.

Therefore, this package provides also a module called HLTL1MuonMatcher that will read the corresponding pat::TriggerObjectStandAlone produced by PAT, but perform the matching using the L1MuonMatcherAlgo. The output will be:

  • A pat trigger association with the matches.
  • A separate collection of fake pat::TriggerObjectStandAlone objects corresponding to the reco object propagated to the second muon station. It has instance label propagatedReco.
  • Another trigger association to the above fake trigger objects (with a null match for reco objects that fail the propagation). It has instance label propagatedReco.
  • If writeExtraInfo is set to cms.bool(True), it will produce two ValueMap<float> with the ΔR and Δφ between the propagated reco and the L1. Instance labels are deltaR, deltaPhi.

The configuration a combination of the parameters of the L1MuonMatcher and the usual parameters from PAT trigger matching modules to HLT path in question and for ambiguity resolution.
Parameters specific to L1 matching are:

name type required notes
src cms.InputTag yes the reco objects, something that can be read through =View
matched cms.InputTag yes the collection of PAT Trigger objects, usually called patTrigger
useTrack cms.string yes one of none, tracker, muon, global. It defines the type of the track that should be extracted from the Candidate object.
none means to just make a fake track with the candidate's momentum and charge.
Note that objects of type RecoChanrgedCandidate , like the ones used by L2 or L3 muon reco, will believe they have a tracker track, not a muon or global one.
useState cms.string sometimes Define which state to start from: atVertex (aka at PCA), innermost, outermost.
Required if useTrack is not none
Note that anything except atVertex will require the reco::TrackExtra to be available.
useSimpleGeometry cms.bool yes If true, the muon station geometry is taken directly from the DetLayers (i.e. one cylinder and two disks),
if false, a further propagation to the individual chambers will be attempted.
The Recommended value is true.
maxDeltaR cms.double yes Maximum allowed deltaR between the L1 candidate and the extrapolated reco object
maxDeltaPhi cms.double no If you want, you can specify also a maximum deltaPhi
sortByDeltaPhi cms.bool no Sort matches by deltaPhi instead of deltaR (if not specified, it will be considered false)
setPropLabel cms.string yes The fake filter label to assign to the pat::TriggerObjectStandAlone made to wrap the 4-momentum of the reco object propagated to the second muon station
Parameters specific to the selection of HLT paths, as for the other PAT Trigger matchers. All parameters are required, but wildcards can be specified in the values.
name type wildcard notes
pathNames cms.vstring '*=, '@' Select only trigger objects used in one of the given trigger paths
filterLabels cms.vstring '*', '@' Select only trigger objects used in one of the given filters
collectionTags cms.vstring '*', '@' Select only trigger objects originating from one of the given collections; note that you have to specify the full tag including the process name (e.g. hltL3MuonCandidates::HLT)
filterIdsEnum cms.vstring '*', '@' Select only trigger objects with one of the given filter IDs assigned to it,
the filter IDs are given as strings according to the enums defined in enum trigger::TriggerObjectType DataFormats/HLTReco/interface/TriggerTypeDefs.h
filterIds cms.vint32 0 Select only trigger objects with one of the given filter IDs assigned to it
andOr cms.bool n/a False means apply the AND of the conditions, True means the OR (it follows the same convention of other HLT modules, strange as might be seems)
resolveAmbiguities cms.bool n/a If True, no more than one reco object can be matched to the same L1 object; precedence is given to the reco ones coming first in the src collection.
Note: the resolveByMatchQuality option is not yet implemented for this module

NOTE: In tag V01-04-00 there was no ambiguity resolution, and selection of the trigger objects was different, controlled by just these two parameters

name type required notes
pathName cms.string yes Name of the HLT path (e.g. HLT_L1MuOpen) to consider
filterLabel cms.string no Name of the filter module; if not specified or empty, it will consider all filters in the path

Code: plugins/HLTL1MuonMatcher.cc (since tag V01-03-00) Config: muonHLTL1Match_cfi.py (since tag V01-03-00)

Associator to trigger objects using HLTDEBUG or HLTMON information

The TriggerMatcherToHLTDebug CMSSW module performs association of reconstructed muons to muon trigger object using the information inside HLTDEBUG or HLTMON. The logic of the matching is:

  1. The track of the reconstructed muon is propagated to station 2
  2. L1 muon are searched for in a ΔR cone around the reco track, and a quality cut is applied to them
  3. For all L1s inside of the cone and passing the cut, the module searches for a L2 seed referencing back to that L1, than a L2 track from that seed, a L3 seed from the L2 track and so on.

The module reports:

  • The number of L1s inside the cone, L1s passing a quality cut, L2 seeds, L2 muons (i.e. full L2 tracks), L2 muons passing filter, L3 seeds, L3 tracker tracks, L3 muons (i.e.g full L3 gobal tracks) and L3 muons passing a filter.
    This information is stored in the event as a set of edm::ValueMap<int> .
  • The references to the L1, L2 and L3 muons. In case of multiple items, the associator will report the reference to the object that goes further in the trigger chain (e.g. if there are two L1s but only the second produces a L2, the L1 reference will be to the second object), or the first object in case they reach the same level (e.g. if there are two L1s both giving a L2 seed and none giving a L2 muon, the first L1 will be reported).
    This information is stored in the event a set of edm::ValueMap<reco::CandidatePtr> .

The module was briefly discussed on hypernews muon-hlt/277.

Code: plugins/TriggerMatcherToHLTDebug.cc (since tag V01-07-00)
Config: triggerMatcherToHLTDebug_cfi.py (since tag V01-07-00)

Configuration file to produce PAT muons with all trigger matches UPDATED

A configuration file is available under MuonAnalysis/MuonAssociators/python/patMuonsWithTrigger_cff.py to produce pat::Muon objects with trigger information, or to add trigger information to existing ones.
Note: before tag V01-10-00 this file was called patMuonsWithTrigger_8E29_cff.

Except for the generic L1 objects, it uses information from the hltTriggerSummaryAOD, so it can see only the trigger filters which report their objects correctly at the end of the HLT

By default, the CFF file produces these matchings:

  • to the nearest L1 object, matched by propagation at the second muon station (ΔR < 0.5)
  • to the nearest L1 passthrough, matched in the same way, and with ambiguity resolution (the same L1 can't be associated to more than one reco muon; the highest pt reco muon has the priority; read below on how to remove this fearure)
  • to the nearest L2 object, matched at the p.c.a. to the beam line (ΔR < 1.2, |ΔpT/pT| < 10; this cuts have not been tuned yet ), with the same ambiguity resolution
  • to the nearest L3 object, matched at the p.c.a. to the beam line (ΔR < 0.5, |ΔpT/pT| < 10; this cuts have not been tuned yet ), with the same ambiguity resolution

Producing the pat::Muons

Producing new pat::Muons with trigger association:
In this simple configuration, pat::Muons will be produced from all reco::Muons inside the muons collection, with all extra features turned off except for the embedding of the tracks. The output collection will be called patMuonsWithTrigger. If you want to customize the producer of pat::Muons (e.g. to add extra information), it can be accessed as process.patMuonsWithoutTrigger.

process.load("MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff")
process.p = cms.Path(process.patMuonsWithTriggerSequence)

Changing the input collection of reco::Muons:

process.load("MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff")
from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff import changeRecoMuonInput
changeRecoMuonInput(process, "myRecoMuons")

Add trigger information to existing pat::Muons:
If you already have your pat::Muons and you want to add trigger matching information to them, this can be easily done

process.load("MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff")
from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff import useExistingPATMuons
useExistingPATMuons(process, "myPatMuons")
process.p = cms.Path(
      process.yourOwnPatMuonSequence *    ## whatever produces the myPatMuons
      process.patMuonsWithTriggerSequence
)
However, there is some subtlety related to the matching with the L1 objects (not the HLT objects corresponding to L1 passthroughs)
  • If you don't need that matching, and you're ok with the matching to L1 passthroughs, the above recipe is ok
  • Wrench, tools otherwise just wait because the current cff is not correct. it will be fixed as soon as possible (but if you need it even sooner, send me a mail)

Changing L1 match cuts to allow single-station triggers in the endcaps:
By default, only muons that reach the second muon station are matched to L1

process.load("MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff")
from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff import useL1MatchingWindowForSinglets
useL1MatchingWindowForSinglets(process)

Adding special di-muon triggers:
StopThis is no longer needed: since tag V01-10-00 all di-muon triggers are also enabled by default Stop
In tags before V01-10-00 the default configuration handles already most of the di-muon triggers, but it does not work for the track leg of the HLT_Mu<X>_Track<X>_Jpsi and HLT_Mu<X>_TkMu<X>_Jpsi triggers. In order to add also that information, add two lines to your config file

process.load("MuonAnalysis.MuonAssociators.patMuonsWithTrigger_8E29_cff")
from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_8E29_cff import addDiMuonTriggers   ## <<===
addDiMuonTriggers(process)   ## <<===
process.p = cms.Path(process.patMuonsWithTriggerSequence)

Removing the ambiguity resolution:
If your signal contains only single muons but you have problems with ghosts, you can switch off the ambiguity resolution

from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff switchOffAmbiguityResolution
switchOffAmbiguityResolution(process)

Changing trigger process name:
In some MC samples, the trigger process is not called HLT. In particular, all 3.5.X samples obtained by a re-digi of 3.1.X samples have REDIGI as trigger process name. You can change this by using the changeTriggerProcessName helper.

from MuonAnalysis.MuonAssociators.patMuonsWithTrigger_cff changeTriggerProcessName
changeTriggerProcessName(process, "REDIGI", oldProcessName="HLT")   # REDIGI for most Spring10 3.5.X MC, REDIGI36X for Summer10 3.6.X MC, REDIGI37X for Summer10 3.7.X MC
If for any bizarre reason you have to call changeTriggerProcessName more than once in the same cfg, each time except the first one you have to add an extra argument oldProcessName="something" specifying the name you have switched to in the (N-1)th attempt...

Reading the HLT information from the pat::Muons

The HLT information can be accessed from the pat::Muon objects by using triggerObjectMatchesByPath("path_name") that will return the 4-vector of the corresponding trigger objects (as pat::TriggerObjectStandAlone), or an empty collection if the muon didn't fire the trigger. The issues with triggerObjectMatchesByPath reporting false positives is fixed in CMSSW_4_1_3 and later.

Alternatively, if you know the internals of the HLT tables, you can use triggerObjectMatchesByFilter and triggerObjectMatchesByCollection to access this information. This is necessary for asymmetric di-muon triggers, as described below.

IDEA! Note about L1 passthrough matching in general:
The matching with L1 passthrough requires for the muon to be successfully propagated to station 2, and this can fail if the muon is too soft; normally, if the muon can't be propagated then we don't expect it to fire any trigger, so this is not a problem. However, if you want, you can check explicitly if the propagation succeeded or not (this can be also useful if you're interested in comparing the L1 (η, φ) with the impact point of the propagated reco. track):

pat::TriggerObjectStandAloneCollection propagated = mu->triggerObjectMatchesByFilter("propagatedToM2");
if (propagated.empty()) {
   std::cout << "Muon didn't reach station 2, according to CMS.SteppingHelixPropagator" << std::end;
} else {
   std::cout << "Propagation succeeeded; eta = " << propagated[0].eta() << ", phi = " << propagated[0].phi() << std::end;
}

HELP Note about "muon plus track" and "muon plus tracker muon" triggers (e.g. MuX_TkMu0_Jpsi):
For these triggers, the triggerObjectMatchesByPath or triggerObjectMatchesByFilter method will return a collection containing both the matches to the L3 muons and the matches to the CKF tracks. In order to see to exactly what trigger was matched, combine this information with the one from triggerObjectMatchesByCollection

bool matchedMu3 = false, matchedTrack = false;
pat::TriggerObjectStandAloneCollection mu0tkMatch = mu->triggerObjectMatchesByFilter("hltMu0TrackJpsiTrackMassFiltered"); // or "hltMu3TkMuJpsiTkMuMassFiltered"
for (unsigned k = 0; k < mu0tkMatch.size(); ++k) {
   if (mu0tkMatch[k].hasCollection("hltL3MuonCandidates") matchedMu = true;
   if (mu0tkMatch[k].hasCollection("hltMuTrackJpsiCtfTrackCands") matchedTrack = true; // or hltMuTkMuJpsiTrackerMuonCands for TkMu0 paths
}
For string-based selectors, e.g. those used in T&P, you can rephrase the above so that it's a single expression. e.g. for matching to the muon and the track part of Mu0_Track0 do
Mu0_Track0_Jpsi_MU = cms.string("!triggerObjectMatchesByCollection('hltL3MuonCandidates').empty() && "+
                                " triggerObjectMatchesByCollection('hltL3MuonCandidates').at(0).hasFilterLabel('hltMu0TrackJpsiTrackMassFiltered')"),
Mu0_Track0_Jpsi_TK = cms.string("!triggerObjectMatchesByCollection('hltMuTrackJpsiCtfTrackCands').empty() && "+
                                " triggerObjectMatchesByCollection('hltMuTrackJpsiCtfTrackCands').at(0).hasFilterLabel('hltMu0TrackJpsiTrackMassFiltered')"),

Note: before CMSSW_3_8_X instead of <obj>.hasCollection("<label>") you have to use <obj>.collection() = "<label>::HLT"

HELP Example macro to read trigger matching:
An example macro that prints out the triggers for a is HeavyFlavorAnalysis/Onia2MuMu/test/printTriggers.cxx. It needs to be compiled in order to check the MuX_Track0_JPsi paths (CINT is not able to understand the code). You can execute it with

root.exe -b -l -q myfile.root printTriggers.cxx+ | tee trigger.txt
Note that the trigger can print the same trigger more than once for a given muon; this is a feature in how trigger matches are embedded but it's not a problem and does not mean that that muon was matched to multiple trigger primitives from the same path.

Reading the L1 information from the pat::Muons

Simulation-related code

Muon Classification By Hits NEW

Matching to the simulation truth at the hit level is possible when working on the GEN-SIM-RECO data tier. The general tools for this studies are in the SimMuon/MCTruth package, and documented in the twikipage SWGuideMuonTrackMCTruth (TODO it misses the new way of runnig on GEN-SIM-RECO without GEN-SIM-RAW without re-digi).

The MuonAssociators package provides a simple module that uses the muon associator by hits tools to classify reconstructed muons in different categories (primary muons, heavy flavours, light flavours and decays, ...). The module is an EDProducer, that writes to the event an edm::ValueMap<int> with the classification code associated to each reco::Muon object; it also writes in other ValueMaps other information that can be used to delve further in the origin of those muons (pdgIds of the mother and grand-mother of the muon, ρ and z at which the muon was produced, ...)

The classification is defined as:

Number Meaning
0 Muon that can't be traced to any simulated hit; noise?
1 Associated to a particle which is not a muon
2 Associated to a muon that comes from a decay in flight or a light flavoured hadron
3 Associated to a prompt muon from the decay of a heavy flavoured hadron (beauty and charm) or a tau
4 Associated to a prompt muon which does not come from hadrons or taus (e.g muons from W/Z/γ*)
The number will be positive for normal muons, and negative for ghosts. Ghost muons are defined as those reconstructed associated to a simulated particle which however matches back to another reconstructed muon; that is, duplicate muons whose hits are already used by a better reconstructed muon (better = more hits shared with the simulated track); note that if you apply some muon selection cuts after the classification, the definition of ghost can become wrong, especially with cuts that have some arbitration inside (e.g. one sim muon gives two tracker muons; you classify before any cuts, the tracker muon 1 has more shared hits with the tracking particle so it's classified as "real" and the muon 2 becomes "ghost"; then you apply arbitration and only muon 2 survives; you've got an event that has only a ghost, which is wrong. the correct solution is to apply the selection before, then you have only one tracker muon and no ghosts)

Code: MuonAnalysis/MuonAssociators/plugins/MuonMCClassifier.cc
Configuration: MuonAnalysis.MuonAssociators.muonClassificationByHits_cfi Example: MuonAnalysis/MuonAssociators/test/MuonMCClassifier/patMuonsWithClassByHits_cfg.py

Tags needed until the code is integrated in release

If you're using a release older than CMSSW_3_8_0_patch2 then you need extra tags in addition to MuonAnalysis/MuonAssociators

For muons reconstructed in CMSSW_3_5_X:

    SimGeneral/TrackingAnalysis V04-01-05
    SimTracker/TrackAssociation V01-08-17
    SimMuon/MCTruth  V02-05-00-03
    MuonAnalysis/MuonAssociators V01-08-08

For muons reconstructed in CMSSW_3_6_X and later:

    SimGeneral/TrackingAnalysis V04-01-05  (already in CMSSW_3_6_1_patch4 or later)
    SimTracker/TrackAssociation V01-08-17
    SimMuon/MCTruth V02-06-01
    MuonAnalysis/MuonAssociators V01-08-08

Example

1) Run the module to perform the association:

  • Load MuonAnalysis.MuonAssociators.muonClassificationByHits_cfi
    • By default, this configuration creates four modules classByHitsTM (all tracker muons), classByHitsTMLSAT (tracker muons last station ang tight), classByHitsGlb (global muons) and classByHitsSta (standAlone muons) that classify muons according to the hits in the matched segments and of the global track (falling back to the standalone track if the global track doesn't have muon hits).
    • If you're using a muon collection other than muons, replace the parameter muons of the four modules with your InputTag
  • Run the muonClassificationByHits sequence before your analyzer

If you want to apply classification to another set of muons with a different selection (e.g. "tracker muosn arbitrated" or "global && tracker") and you want to be sure that you get the ghosts defined correctly you should make a clone of one of those four modules, change the muonPreselection parameter to your selection cut and run it after the muonClassificationByHits sequence. See for example patMuonsWithClassByHits_cfg.py.
This is important for selections that can have some effect on arbitration, because the only situation in which applying the selection afterwards breaks the classification is when you have a muon reconstructed multiple times and your selection cut could kill the one classified as "real" and let the "ghost" survive.

2A) Use the classification by hits in your analyzer code

You can use the classification code easily in your EDAnalyzer. The only thing you have to take care is that if you're selecting muons by value (e.g. using MuonSelector instead of MuonRefSelector) you have to run the classifier on the already selected muons.

Handle<View<reco::Muon> > muons; iEvent.getByLabel(src_, muons);
Handle<ValueMap<int> > classif; iEvent.getByLabel(classif_, classif);

for (size_t i = 0, n = muons->size(); i < n; ++i) {
   edm::RefToBase<reco::Muon> muRef = muons->refAt(i);
   int origin = (*classif)[muRef];
   if (origin < 0) {
       std::cout << "This is a ghost!" << std::endl;
   } else {
      switch (origin) {   
          case 0: std::cout << "Unmatched (??)" << std::endl; break;
          case 1: std::cout << "Fake" << std::endl; break;
          case 2: std::cout << "Light flavour or decay" << std::endl; break;          
          case 3: std::cout << "Heavy flavour or tau" << std::endl; break;
          case 4: std::cout << "Primary muon" << std::endl; break;
      }
   }
}

2B) Write the classification by hits into pat::Muons
The information of the classification can be easily embedded inside pat::Muons to be accessible e.g. in FWLite. If you PATMuonProducer is called patMuons, then you can just put in your python config file

from MuonAnalysis.MuonAssociators.muonClassificationByHits_cfi import addUserData as addClassByHits
addClassByHits(process.patMuons)
Afterwards, you'll be able to access the classification as muon.userInt("classByHitsTM") for the classification done using tracker muon segments, or muon.userInt("classByHitsGlb") for the classification done using the hits inside of the muon track.

Further information TODO

The module produces also other ValueMaps with more detailed information. It will be documented here ASAP, but in the meantime you can refer to the code.

GenParticle Associator by Helix Pulls

Matches a candidate containing a track to some Candidate or GenParticle using the pulls of the track parameters at the point of closest approach. A preselection by deltaR is used to reduce the combinatorics.

The user can choose if off-diagonal terms of the covariance matrix should be used or not, and if the two vertex-related parameters dxy and dsz should be included in the pull or not.

See also this hypernews thread

Algorithm

The algorithm object takes care of the matching itself, and can handle also multiple matches.

The two high level match methods provided to the user are:

  • pair<int,float> match(RecoCandidate mu, vector<GenParticle> gens, vector<uint8_t> good)
    Try to match the muon to the GenParticles in gens, considering only those for which good[i] is true.
    On a success, it will return the index and the pull; on a failure, it will return the index and the pull. On a failure, it will return (-1,9e9).
  • void matchMany(RecoCandidate mu, vector<GenParticle> gens, vector<uint8_t> good, vector<pair<int,float> > & toFill)
    Try to match the muon to the GenParticles in gens, considering only those for which good[i] is true, and allowing multiple matches.
    Matches sorted by pulls will be pushed back into the vector toFill passed by reference; the toFill vector must be empty at the beginning.
two similar methods are provided taking as input bare reco::Track objects
  • pair<int,float> match(reco::Track tk, vector<GenParticle> gens, vector<uint8_t> good)
  • void matchMany(reco::Track tk, vector<GenParticle> gens, vector<uint8_t> good, vector<pair<int,float> > & toFill)

There are also a few lower level methods, used by the algorithm itself but also declared public as they can be useful to more advanced users.

  • pair<bool,float> match(Track tk, Candidate mc, AlgebraicSymMatrix55 invertedCovariance)
    Match the track helix to the mc particle; it will use the inverse of the helix covariance matrix passed as parameter to compute the pull. Return is (true,pull) on a success, (false,9e9) on a failure.
  • void fillInvCov(Track tk, AlgebraicSymMatrix55 & invertedCovariance)
    Extract the covariance matrix from the track object, possibly take only the diagonal part and/or remove the vertex variables, and invert the matrix.

Example usage:

edm::Handle<std::vector<reco::GenParticle> > gens;
   iEvent.getByLabel(gens_, gens);
   std::vector<uint8_t> good; good.reserve(gens->size());
   for (std::vector<reco::GenParticle>::const_iterator it = gens->begin(), ed = gens->end(); it != ed; ++it) { 
       good.push_back( (it->pt() > 1)  && (it->status() == 1) );
   }

   const reco::Muon &amp; mu = ...
   std::pair<int,float> match = matcher_.match(mu, *gens, good);
   if (match.first != -1) {
      std::cout << "Matched with GenParticle[" << match.first << "], pull is " << match.second << std::endl;
   }

The matcher takes the following configuration parameters:

name type *notes
track cms.string Which reco::RecoCandidate method to use to get the track: one of tracker, standAloneMuon, combinedMuon.
Note that objects of type RecoChanrgedCandidate , like the ones used by L2 or L3 muon reco, will believe they have a tracker track, not a muon or global one.
maxDeltaR cms.double A preselection cut on the deltaR to reduce the number of pairs to be considered
maxPull cms.double The cut to be applied on the pull
useVertexVariables cms.bool Include dxy, dsz in the pull computation.
diagonalElementsOnly cms.bool Ignore off-diagonal elements of the covariance matrix when computing the pulls.
The diagonal terms are q/p, theta, phi, dxy, dsz.

Code: interface/MatcherByPullsAlgorithm.h, src/MatcherByPullsAlgorithm.cc (since tag V00-00-01)

Warning, important To be done:

  • Find out sensible values for the cut on the pull, and check the difference in different modes (with/without vertex, using also off-diagonal terms)
  • Possibly apply some propagation to PCA for particles produced at a radius larger than few centimeters

CMSSW Module (for Candidate objects)

The module runs the matcher and produces and Associations to GenParticles and a ValueMaps with the pulls. This data can be used from an EDAnalyzer, or copied into a pat::Muon, or pat::GenericParticle.

In addition to the parameters from the algorithm, it requires also

name type *notes
src cms.InputTag The collection of the input objects to match, which must inherit from reco::RecoCandidates (e.g. reco::Muons)
matched cms.InputTag The collection of the GenParticles to match against. It must be a reco::GenParticleCollection (usually it's genParticles)
matchedSelector cms.string A generic cut on the GenParticles

Example configuration:

process.matchByPulls = cms.EDProducer("MatcherByPulls",
    src = cms.InputTag("muons"),
    matched         = cms.InputTag("genParticles"),
    matchedSelector = cms.string("status == 1"),
    track = cms.string("track"),
    maxDeltaR = cms.double(0.3),
    maxPull   = cms.double(200),
    diagonalElementsOnly = cms.bool(True),
    useVertexVariables   = cms.bool(True),
)

Example using pat::Muons
To use the matcher for a pat::Muon

process.patDefaultSequence.replace(process.muonMatch, process.matchByPulls)
process.allLayer1Muons.userData.userFloats.src += [ cms.InputTag("matchByPulls","pulls") ]
process.allLayer1Muons.genParticleMatch = cms.InputTag("matchByPulls") 
To read back the information in your analyzer code, you can then check if the muon matches simply by
reco::GenParticleRef gen = mu.genParticleRef();
if (gen.isNull()) { 
   // sorry, no match 
} else { 
   float genPt = gen->pt(); 
   float pull    = mu->userFloat("matcherByPulls:pulls");
}

Example of standalone usage

// Get the matches
Handle<Association<reco::GenParticleCollection> > matches;
iEvent.getByLabel(matcherTag_, matches);

// Get the pulls
Handle<ValueMap<float> > pulls;
iEvent.getByLabel(edm::InputTag(matcherTag_.label(), "pulls", matcherTag_.process()), pulls);

// Get the muon
Handle<View<reco::Muon> > muons;
iEvent.getByLabel(muons_, muons);

for (size_t i = 0, n = muons->size(); i < n; ++i) {
    Ptr<reco::Muon> muPtr = muons->ptrAt(i);       // RefToBase would work too, using 'refAt(i)'
    reco::GenParticleRef gen = (*matches)[muPtr];
    if (gen.isNonnull()) {
        float genPt = gen->pt(); 
        float pull    = (*pulls)[muPtr];
    }
}

Code: plugins/MatcherByPulls.cc (since tag V00-00-01)

CMSSW Module (for bare reco::Track objects)

Since tag V01-02-00 there is also an EDProducer TrackMatcherByPulls that can take as input a collection of reco::Track objects instead of the Candidates. It takes the same input parameters as the MatcherByPulls, and produces an output in the same format.

Example configuration

process.matchByPulls = cms.EDProducer("TrackMatcherByPulls",
    src = cms.InputTag("generalTracks"),
    matched         = cms.InputTag("genParticles"),
    matchedSelector = cms.string("status == 1"),
    track = cms.string("track"),
    maxDeltaR = cms.double(0.3),
    maxPull   = cms.double(200),
    diagonalElementsOnly = cms.bool(True),
    useVertexVariables   = cms.bool(True),
)

Example code

// Get the matches
edm::Handle<edm::Association<reco::GenParticleCollection> > matches;
iEvent.getByLabel(matcherTag_, matches);

// Get the pulls
edm::Handle<edm::ValueMap<float> > pulls;
iEvent.getByLabel(edm::InputTag(matcherTag_.label(), "pulls", matcherTag_.process()), pulls);

// Get the muon
edm::Handle<edm::View<reco::Track> > tracks;
iEvent.getByLabel(tracks_, tracks);

for (size_t i = 0, n = tracks->size(); i < n; ++i) {
    edm::RefToBase<reco::Track> tk = tracks->refAt(i);   // Ptr would work too, using 'ptrAt(i)'
    reco::GenParticleRef gen = (*matches)[tk];
    if (gen.isNonnull()) {
        float genPt = gen->pt(); 
        float pull    = (*pulls)[muPtr];
    }
}

Code related to reconstruction, or generic

Matcher for objects containing tracks

MatcherUsingTracksAlgorithm is a generic tool that can be used to match RECO objects to other RECO objects by looking at the tracks inside them. In particular, it can:

  • Check if two objects refer to the same track
  • Compare the tracks in the two objects at any given state (at vertex, innermost, outermost)
  • Extrapolate then track or momentum of one object to the surface defined by the innermost or outermost state of the track of the other object, and match there.
In the future, it might also be able to perform matching by hits

The usecases are tag and probe and matching, e.g. of offline reconstruction to the HLT one or of tracker tracks to standalone muon tracks.

Algorithm

The algorithm object takes care of the matching, either between two objects or between one object and a list of other objects.

The relevant methods are:

  • bool match(const reco::Candidate &src, const reco::Candidate &matched, float &deltaR, float &deltaLocalPos, float &deltaPtRel, float & chi2) const
  • int match(const reco::Candidate &src, const edm::View<reco::Candidate> & matcheds; float &deltaR, float &deltaLocalPos, float &deltaPtRel, float & chi2) const
Both will try to match the object src against the object or objects in matched; the first method will return true if the two objects match, the second will return the index of the best match or -1 if no match is found.
The values of deltaR, deltaLocalPos, deltaPtRel and chi2 should be set to some large number when calling the method, and are filled with the results of the best match if a match is found.

Auxiliary methods are:

  • init(const edm::EventSetup &) to be called at the beginning of an event or run so that the matcher can retrieve from the EventSetup the necessary services (geometry, magnetic field and propagator)
  • bool hasMetrics() will return true if the algorithm is computing deltaR, deltaLocalPos and deltaPtRel values, or false otherwise
  • bool hasChi2() will return true if the algorithm is computing the chi2

The matcher object can be built from a ParameterSet, with the following configurables:

name type notes
algorithm cms.string How to perform the matching:
byTrackRef equality between references
byDirectComparison compare the track parameters directly (usually the PCA)
byPropagatingSrc, byPropagatingMatched extrapolate one track to the surface defined by the other track.
srcTrack
matchedTrack
cms.string choice of track for candidates src and matched. Values are tracker, muon, global, none (i.e. use the Candidate's charge and momentum)
requireSameCharge cms.bool Require the two tracks to thave the the same charge. This parameter is optional (default is false), and only available since V01-12-01
Additional parameters are needed for all matchers except byTrackRef
name type notes
srcState
matchedState
cms.string Which state to use: atVertex (i.e. at the PCA to the beam line), innermost, outermost
maxDeltaR
(maxDeltaEta,
maxDeltaPhi )
cms.double Cut on the deltaR between the momentum of the two tracks (at whatever state the matching is performed).
Additional cuts on deltaEta and deltaPhi can be optionally specified.
maxDeltaLocalPos cms.double Cut on the distance in between the two positions
maxDeltaPtRel cms.double Cut on (ptsrc - ptmatched>)/(ptmatched)
sortBy cms.string Choose which variable will be used to sort matches and pick the best one in case

Since tag V00-01-01, the module can also compute the compatibility chi2 between two track states This can be activated by setting computeChi2 to cms.bool(True) ; if the parameter is not present, it defaults to false. When the chi2 computation is turned on, other parameters are necessary

name type notes
maxChi2 cms.double Cut on the compatibility chi2 in the matching
chi2DiagonalOnly cms.bool Ignore off-diagonal terms in the covariance matrix when computing the chi2
chi2UsePosition cms.bool Usable only when algorithm is byPropagatingSrc or byPropagatingMatched, includes the local position difference in the chi2; if set to false, only the momentum is used.
chi2UseVertex cms.bool Usable only when algorithm is byDirectComparison , includes the dxy and dsz variables in the chi2; if set to false, only the momentum is used.
chi2MomentumForDxy cms.string Usable only when algorithm is byDirectComparison , determines which of the two track momenta is to be used in the dxy computation; allowed values are src and matched.

Code: interface/MatcherUsingTracksAlgorithm.h, src/MatcherUsingTracksAlgorithm.cc (since tag V00-00-02)

CMSSW Module

The cmssw module, called MatcherUsingTracks, will perform the match of items in one collection of Candidates to items in another one.

The parameters of the module are:

name type notes
src cms.InputTag collections of Candidates objects you want to match (e.g. offline muons, tracker tracks)
matched cms.InputTag collections of Candidates against which src should be matched (e.g. L2 muons or STA tracks)
dontFailOnMissingInput cms.bool If set to true, the module treat a missing collection as if it was empty instead of throwing the usual ProductNotFound exception. This is needed e.g. in case the collection is the output of an HLT module, that might not have been executed in a given event if a filter before it failed.
The parameter is tracked but optional, the default value is false.
writeExtraPATOutput cms.bool If set to true, the will write the match output also as edm::ValueMap<edm::Ptr<pat::CMS.UserData> > (which might require dictionaries not available in the release).
The parameter is tracked but optional, the default value is false.

In addition, it takes all the parameters of the MatcherUsingTracksAlgorithm

The module produces as output a ValueMap<reco::CandidatePtr> to the matched objects, plus some other ValueMaps of ints and floats to hold extra information about the matching:

  • for the matching by track ref, it will produce a ValueMap<int> with instance label matched, containing 0 for unmatched objects, 1 for matched ones.
  • for the other matchings it will produce three ValueMaps of floats, with instance labels deltaR, deltaEta, deltaPhi, deltaLocalPos and deltaPtRel, holding the values for these matches (or 999 in case of no match). If the matcher also computes the chi2, it will also produce a ValueMap<float> for it.

Code: plugins/MatcherUsingTracks.cc (since tag V00-00-02)

Example configurations:

  • To match some particles with tracker tracks (e.g. RecoChargedCandidates or pat::!GenericParticles made from generalTracks) to standalone muons, propagating the tracker track to the innermost muon state (requires the TrackExtra for the muon track, not the tracker one, so it works on AOD)
process.staToTkMatch = cms.EDProducer("MatcherUsingTracks",
    src     = cms.InputTag("selectedLayer1TrackCands"),  # these were made from generalTracks
    matched = cms.InputTag("standAlonePATMuons"),     # subset of PAT muons with a standalone muon track
    algorithm = cms.string("byPropagatingSrc"),
    srcTrack = cms.string("tracker"),
    srcState = cms.string("atVertex"),
    matchedTrack = cms.string("muon"),
    matchedState = cms.string("innermost"),
    maxDeltaR        = cms.double(0.3),  # this makes almost sense with the current CRAFT-based alignment
    maxDeltaLocalPos = cms.double(100), # set it to 1m, so it's not really cutting anything
    maxDeltaPtRel    = cms.double(10),     # same as above, allow for very large mismeasurement in pT
    sortBy           = cms.string("deltaR"),
)
  • Reverse of the above: match standalone muons to tracker tracks, but still do the tracker to muon propagation instead of the muon to tracker one (so there are no issues with muon resolution or B field uncertainties)
process.staToTkMatch = cms.EDProducer("MatcherUsingTracks",
    src     = cms.InputTag("standAlonePATMuons"),
    matched = cms.InputTag("selectedLayer1TrackCands"),
    algorithm = cms.string("byPropagatingMatched"),
    srcTrack = cms.string("muon"),
    srcState = cms.string("innermost"),
    matchedTrack = cms.string("tracker"),
    matchedState = cms.string("atVertex"),
    maxDeltaR        = cms.double(0.3),
    maxDeltaLocalPos = cms.double(100),
    maxDeltaPtRel    = cms.double(10),
    sortBy           = cms.string("deltaR"),
)
  • To match some particles with tracker tracks (e.g. RecoChargedCandidates or pat::!GenericParticles made from generalTracks) to reco muons requiring that they use the very same track.
process.staToTkMatch = cms.EDProducer("MatcherUsingTracks",
    src     = cms.InputTag("allTrackCands"),  # these were made from generalTracks
    matched = cms.InputTag("muons"),     
    algorithm = cms.string("byTrackRef"),
    srcTrack = cms.string("tracker"),
    matchedTrack = cms.string("tracker"),
)
  • Match offline muons to L2 muons, propagating the offline tracker track to the L2 muon. Requires RAW+HLTMON data tier.
process.matchL2Tk = cms.EDProducer("MatcherUsingTracks",
    src     = cms.InputTag("muons"),
    matched = cms.InputTag("hltL2MuonCandidates","","HLT"),
    dontFailOnMissingInput = cms.bool(True), # If HLT doesn't find candidates it won't even produce an empty product!
    algorithm    = cms.string("byPropagatingSrc"),
    srcTrack     = cms.string("tracker"),
    srcState     = cms.string("outermost"),
    matchedTrack = cms.string("tracker"),  ## L2 muons are RecoChargedCandidate, so they think they have a tracker track
    matchedState = cms.string("innermost"),
    maxDeltaLocalPos = cms.double(200),
    maxDeltaR        = cms.double(0.3),
    maxDeltaPtRel    = cms.double(999.0),
    sortBy           = cms.string("deltaR")
)
  • Match offline muons to L3 muons at vertex
process.matchL3Tk = cms.EDProducer("MatcherUsingTracks",
    src     = cms.InputTag("muons"),
    matched = cms.InputTag("hltL3MuonCandidates","","HLT"),
    dontFailOnMissingInput = cms.bool(True), # If HLT doesn't find candidates it won't even produce an empty product!
    algorithm    = cms.string("byDirectComparison"),
    srcTrack     = cms.string("tracker"),
    srcState     = cms.string("atVertex"),
    matchedTrack = cms.string("tracker"),   # same as above, a L3 object thinks it has a tracker track
    matchedState = cms.string("atVertex"),
    maxDeltaLocalPos = cms.double(200),
    maxDeltaR        = cms.double(0.3),
    maxDeltaPtRel    = cms.double(999.0),
    sortBy           = cms.string("deltaR")
)

Some other examples config files are available under UserCode/GPetrucc/test/TrackMatcher; there is also an fwlite macro about how to access some of the matching information UserCode/GPetrucc/test/TrackMatcher/scan.cc

Review status

Reviewer/Editor and Date (copy from screen) Comments
GiovanniPetrucciani - 18 May 2009 L1 matcher
GiovanniPetrucciani - 29 May 2009 Matcher by pulls
GiovanniPetrucciani - 24 Jun 2009 Matcher using tracks (still lacking examples)

Responsible: GiovanniPetrucciani
Last reviewed by: GiovanniPetrucciani

Edit | Attach | Watch | Print version | History: r33 < r32 < r31 < r30 < r29 | Backlinks | Raw View | WYSIWYG | More topic actions
Topic revision: r33 - 2011-04-24 - GiovanniPetrucciani
 
    • Cern Search Icon Cern Search
    • TWiki Search Icon TWiki Search
    • Google Search Icon Google Search

    CMSPublic All webs login

This site is powered by the TWiki collaboration platform Powered by PerlCopyright &© 2008-2023 by the contributing authors. All material on this collaboration platform is the property of the contributing authors.
or Ideas, requests, problems regarding TWiki? use Discourse or Send feedback