Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

make V0EventSelector produce skimmed V0 collections based on their mass #46005

Merged
merged 1 commit into from
Sep 16, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
28 changes: 25 additions & 3 deletions DQM/TrackingMonitorSource/plugins/V0EventSelector.cc
Original file line number Diff line number Diff line change
Expand Up @@ -2,11 +2,13 @@
#include <vector>
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/Common/interface/Ref.h"
#include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
#include "FWCore/Framework/interface/stream/EDFilter.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"

class V0EventSelector : public edm::stream::EDFilter<> {
public:
Expand All @@ -19,25 +21,45 @@ class V0EventSelector : public edm::stream::EDFilter<> {
private:
const edm::EDGetTokenT<reco::VertexCompositeCandidateCollection> vccToken_;
const unsigned int minNumCandidates_;
const double massMin_;
const double massMax_;
};

V0EventSelector::V0EventSelector(const edm::ParameterSet& iConfig)
: vccToken_{consumes<reco::VertexCompositeCandidateCollection>(
iConfig.getParameter<edm::InputTag>("vertexCompositeCandidates"))},
minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")} {}
minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")},
massMin_{iConfig.getParameter<double>("massMin")},
massMax_{iConfig.getParameter<double>("massMax")} {
produces<reco::VertexCompositeCandidateCollection>();
}

bool V0EventSelector::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
edm::Handle<reco::VertexCompositeCandidateCollection> vccHandle;
iEvent.getByToken(vccToken_, vccHandle);

return vccHandle->size() >= minNumCandidates_;
auto filteredVCC = std::make_unique<reco::VertexCompositeCandidateCollection>();

for (const auto& vcc : *vccHandle) {
if (vcc.mass() >= massMin_ && vcc.mass() <= massMax_) {
filteredVCC->push_back(vcc);
}
}

bool pass = filteredVCC->size() >= minNumCandidates_;
iEvent.put(std::move(filteredVCC));

return pass;
}

void V0EventSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("vertexCompositeCandidates", edm::InputTag("generalV0Candidates:Kshort"));
desc.add<unsigned int>("minCandidates", 1); // Change '1' to your desired minimum number of candidates
desc.add<unsigned int>("minCandidates", 1)->setComment("Minimum number of candidates required");
desc.add<double>("massMin", 0.0)->setComment("Minimum mass in GeV");
desc.add<double>("massMax", std::numeric_limits<double>::max())->setComment("Maximum mass in GeV");
descriptions.addWithDefaultLabel(desc);
}

// Define this module as a plug-in
DEFINE_FWK_MODULE(V0EventSelector);
Original file line number Diff line number Diff line change
Expand Up @@ -94,7 +94,7 @@
ttbarEventSelector = cms.EDFilter("ttbarEventSelector")
ttbarTracks = cms.EDProducer("TtbarTrackProducer")

# Added modules for V0Monitoring
# Added modules for V0Monitoring (standard)
KshortMonitor = v0Monitor.clone()
KshortMonitor.FolderName = "StandaloneTrackMonitor/V0Monitoring/Ks"
KshortMonitor.v0 = "generalV0Candidates:Kshort"
Expand All @@ -108,6 +108,18 @@
LambdaMonitor.histoPSet.massPSet = cms.PSet(nbins = cms.int32(100),
xmin = cms.double(1.050),
xmax = cms.double(1.250))

# Added modules for V0Monitoring (for restricted mass candidates)
SelectedKshortMonitor = KshortMonitor.clone(
FolderName = "StandaloneTrackMonitor/V0Monitoring/SelectedKs",
v0 = "KShortEventSelector"
)

SelectedLambdaMonitor = LambdaMonitor.clone(
FolderName = "StandaloneTrackMonitor/V0Monitoring/SelectedLambda",
v0 = "LambdaEventSelector"
)

##################
# For MinBias
##################
Expand Down Expand Up @@ -171,31 +183,31 @@
* KShortEventSelector
* KshortTracks
* standaloneTrackMonitorK0
* KshortMonitor)
* SelectedKshortMonitor)

standaloneValidationK0sMC = cms.Sequence(
hltPathFilter
* selectedPrimaryVertices
* KShortEventSelector
* KshortTracks
* standaloneTrackMonitorK0MC
* KshortMonitor)
* SelectedKshortMonitor)

standaloneValidationLambdas = cms.Sequence(
hltPathFilter
* selectedPrimaryVertices
* LambdaEventSelector
* LambdaTracks
* standaloneTrackMonitorLambda
* LambdaMonitor)
* SelectedLambdaMonitor)

standaloneValidationLambdasMC = cms.Sequence(
hltPathFilter
* selectedPrimaryVertices
* LambdaEventSelector
* LambdaTracks
* standaloneTrackMonitorLambdaMC
* LambdaMonitor)
* SelectedLambdaMonitor)

##################
# For ZtoEE
Expand Down
17 changes: 13 additions & 4 deletions DQM/TrackingMonitorSource/python/V0Selections_cfi.py
Original file line number Diff line number Diff line change
@@ -1,12 +1,21 @@
from DQM.TrackingMonitorSource.v0EventSelector_cfi import *
from DQM.TrackingMonitorSource.v0VertexTrackProducer_cfi import *

KShortEventSelector = v0EventSelector.clone()
KShortEventSelector = v0EventSelector.clone(
vertexCompositeCandidates = "generalV0Candidates:Kshort",
massMin = 0.47,
massMax = 0.52,
)

LambdaEventSelector = v0EventSelector.clone(
vertexCompositeCandidates = "generalV0Candidates:Lambda"
vertexCompositeCandidates = "generalV0Candidates:Lambda",
massMin = 1.11,
massMax = 1.128
)

KshortTracks = v0VertexTrackProducer.clone()
KshortTracks = v0VertexTrackProducer.clone(
vertexCompositeCandidates = "KShortEventSelector"
)
LambdaTracks = v0VertexTrackProducer.clone(
vertexCompositeCandidates = "generalV0Candidates:Lambda"
vertexCompositeCandidates = "LambdaEventSelector"
)