b Tag Analysis with a CMSSW EDAnalyzer in CMSSW 12.0.X
Most physics analysis will use this method, as it is the most general and flexible way to access and save event information. One should start by learning how to create a generic EDAnalyzer and configuration files, as described
here in the workbook.
In summary, to access b-tag information in each event, one can edit your EDAnalyzer code in the following ways. Add the following header file at the beginning of the code:
#include "DataFormats/BTauReco/interface/JetTag.h"
Add the following line below the include statements:
using reco::JetTagCollection;
In the EDAnalyzer class definition section, define a token to identify the b-tagging discriminator you want to study:
edm::EDGetTokenT<JetTagCollection> bTagToken_;
Finally, we access the b-tag information by adding the following lines in the
EDAnalyzer::analyze()
portion of the code:
// Get b-tag information
edm::Handle<JetTagCollection> bTagHandle = iEvent.getHandle(bTagToken_);
const reco::JetTagCollection &bTags = *(bTagHandle.product());
// Loop over jets and study b-tag info
for(unsigned int i = 0; i != bTags.size(); i++){
printf("[bTagInfo] (pt, eta, phi) jet = %7.2f %6.3f %6.3f |
b-tag descriminator = %6.3f \n",
bTags[i].first->pt(), bTags[i].first->eta(),
bTags[i].first->phi(), bTags[i].second() ) ;
}
Next, we need to add the name of the b-tag input tag to the configuration file. This is done through the variable
bTags
(as defined in the above portion of the code). Add the following line to the configuration code:
bTags = cms.untracked.InputTag('myChosenBTag')
where
myChosenBTag
can be the module name of any b-tag discriminator one wants to analyze. Some examples include
pfCombinedSecondaryVertexV2BJetTags
,
pfJetProbabilityBJetTags
,
softPFMuonBJetTags
, and more. You can find all of the possible tags in a given data file by exploring the output of the command
edmDumpEventContent /path/to/mySourceDataFile.root
Next, we need to edit the
Build.xml
file associated with the EDAnalyzer. Add the line
<user name=DataFormats/BTauReco>
The output of your EDAnalyzer should look like this :
04-Feb-2022 22:26:01 CET Initiating request to open file file:/afs/cern.ch/cms/Tutorials/workbook_twiki2021/MinBias_pythia8_14TeV_100events.root
04-Feb-2022 22:26:07 CET Successfully opened file file:/afs/cern.ch/cms/Tutorials/workbook_twiki2021/MinBias_pythia8_14TeV_100events.root
Begin processing the 1st record. Run 1, Event 1, LumiSection 1 on stream 0 at 04-Feb-2022 22:26:22.670 CET
...
Begin processing the 34th record. Run 1, Event 32, LumiSection 1 on stream 0 at 04-Feb-2022 22:26:24.932 CET
analyzing event run: 1 lumi: 1 event: 32
[bTagInfo] (pt, eta, phi) jet = 51.85 0.323 2.104 | b-tag descriminator = 0.130
[bTagInfo] (pt, eta, phi) jet = 43.04 -0.775 -1.651 | b-tag descriminator = 0.080
[bTagInfo] (pt, eta, phi) jet = 26.93 0.123 -1.657 | b-tag descriminator = 0.250
[bTagInfo] (pt, eta, phi) jet = 21.11 -1.077 1.221 | b-tag descriminator = 0.361
...
Begin processing the 99th record. Run 1, Event 99, LumiSection 1 on stream 0 at 04-Feb-2022 22:26:26.543 CET
analyzing event run: 1 lumi: 1 event: 99
Begin processing the 100th record. Run 1, Event 95, LumiSection 1 on stream 0 at 04-Feb-2022 22:26:26.548 CET
analyzing event run: 1 lumi: 1 event: 95
[bTagInfo] (pt, eta, phi) jet = 23.35 -0.934 -0.067 | b-tag descriminator = 0.362
[bTagInfo] (pt, eta, phi) jet = 21.11 -0.585 -0.919 | b-tag descriminator = 0.144
--
ClaytonBennett - 2022-02-21