Skip to content

Commit c891039

Browse files
authored
Merge pull request #46005 from mmusich/mm_dev_V0EventSelector
make `V0EventSelector` produce skimmed V0 collections based on their mass
2 parents e315d8d + 43db8e1 commit c891039

File tree

3 files changed

+55
-12
lines changed

3 files changed

+55
-12
lines changed

DQM/TrackingMonitorSource/plugins/V0EventSelector.cc

Lines changed: 25 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -2,11 +2,13 @@
22
#include <vector>
33
#include "DataFormats/Common/interface/Handle.h"
44
#include "DataFormats/Common/interface/View.h"
5+
#include "DataFormats/Common/interface/Ref.h"
56
#include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
67
#include "FWCore/Framework/interface/stream/EDFilter.h"
78
#include "FWCore/Framework/interface/Event.h"
89
#include "FWCore/Framework/interface/MakerMacros.h"
910
#include "FWCore/ParameterSet/interface/ParameterSet.h"
11+
#include "FWCore/Utilities/interface/InputTag.h"
1012

1113
class V0EventSelector : public edm::stream::EDFilter<> {
1214
public:
@@ -19,25 +21,45 @@ class V0EventSelector : public edm::stream::EDFilter<> {
1921
private:
2022
const edm::EDGetTokenT<reco::VertexCompositeCandidateCollection> vccToken_;
2123
const unsigned int minNumCandidates_;
24+
const double massMin_;
25+
const double massMax_;
2226
};
2327

2428
V0EventSelector::V0EventSelector(const edm::ParameterSet& iConfig)
2529
: vccToken_{consumes<reco::VertexCompositeCandidateCollection>(
2630
iConfig.getParameter<edm::InputTag>("vertexCompositeCandidates"))},
27-
minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")} {}
31+
minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")},
32+
massMin_{iConfig.getParameter<double>("massMin")},
33+
massMax_{iConfig.getParameter<double>("massMax")} {
34+
produces<reco::VertexCompositeCandidateCollection>();
35+
}
2836

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

33-
return vccHandle->size() >= minNumCandidates_;
41+
auto filteredVCC = std::make_unique<reco::VertexCompositeCandidateCollection>();
42+
43+
for (const auto& vcc : *vccHandle) {
44+
if (vcc.mass() >= massMin_ && vcc.mass() <= massMax_) {
45+
filteredVCC->push_back(vcc);
46+
}
47+
}
48+
49+
bool pass = filteredVCC->size() >= minNumCandidates_;
50+
iEvent.put(std::move(filteredVCC));
51+
52+
return pass;
3453
}
3554

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

64+
// Define this module as a plug-in
4365
DEFINE_FWK_MODULE(V0EventSelector);

DQM/TrackingMonitorSource/python/TrackingDataMCValidation_Standalone_cff.py

Lines changed: 17 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -94,7 +94,7 @@
9494
ttbarEventSelector = cms.EDFilter("ttbarEventSelector")
9595
ttbarTracks = cms.EDProducer("TtbarTrackProducer")
9696

97-
# Added modules for V0Monitoring
97+
# Added modules for V0Monitoring (standard)
9898
KshortMonitor = v0Monitor.clone()
9999
KshortMonitor.FolderName = "StandaloneTrackMonitor/V0Monitoring/Ks"
100100
KshortMonitor.v0 = "generalV0Candidates:Kshort"
@@ -108,6 +108,18 @@
108108
LambdaMonitor.histoPSet.massPSet = cms.PSet(nbins = cms.int32(100),
109109
xmin = cms.double(1.050),
110110
xmax = cms.double(1.250))
111+
112+
# Added modules for V0Monitoring (for restricted mass candidates)
113+
SelectedKshortMonitor = KshortMonitor.clone(
114+
FolderName = "StandaloneTrackMonitor/V0Monitoring/SelectedKs",
115+
v0 = "KShortEventSelector"
116+
)
117+
118+
SelectedLambdaMonitor = LambdaMonitor.clone(
119+
FolderName = "StandaloneTrackMonitor/V0Monitoring/SelectedLambda",
120+
v0 = "LambdaEventSelector"
121+
)
122+
111123
##################
112124
# For MinBias
113125
##################
@@ -171,31 +183,31 @@
171183
* KShortEventSelector
172184
* KshortTracks
173185
* standaloneTrackMonitorK0
174-
* KshortMonitor)
186+
* SelectedKshortMonitor)
175187

176188
standaloneValidationK0sMC = cms.Sequence(
177189
hltPathFilter
178190
* selectedPrimaryVertices
179191
* KShortEventSelector
180192
* KshortTracks
181193
* standaloneTrackMonitorK0MC
182-
* KshortMonitor)
194+
* SelectedKshortMonitor)
183195

184196
standaloneValidationLambdas = cms.Sequence(
185197
hltPathFilter
186198
* selectedPrimaryVertices
187199
* LambdaEventSelector
188200
* LambdaTracks
189201
* standaloneTrackMonitorLambda
190-
* LambdaMonitor)
202+
* SelectedLambdaMonitor)
191203

192204
standaloneValidationLambdasMC = cms.Sequence(
193205
hltPathFilter
194206
* selectedPrimaryVertices
195207
* LambdaEventSelector
196208
* LambdaTracks
197209
* standaloneTrackMonitorLambdaMC
198-
* LambdaMonitor)
210+
* SelectedLambdaMonitor)
199211

200212
##################
201213
# For ZtoEE
Lines changed: 13 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -1,12 +1,21 @@
11
from DQM.TrackingMonitorSource.v0EventSelector_cfi import *
22
from DQM.TrackingMonitorSource.v0VertexTrackProducer_cfi import *
33

4-
KShortEventSelector = v0EventSelector.clone()
4+
KShortEventSelector = v0EventSelector.clone(
5+
vertexCompositeCandidates = "generalV0Candidates:Kshort",
6+
massMin = 0.47,
7+
massMax = 0.52,
8+
)
9+
510
LambdaEventSelector = v0EventSelector.clone(
6-
vertexCompositeCandidates = "generalV0Candidates:Lambda"
11+
vertexCompositeCandidates = "generalV0Candidates:Lambda",
12+
massMin = 1.11,
13+
massMax = 1.128
714
)
815

9-
KshortTracks = v0VertexTrackProducer.clone()
16+
KshortTracks = v0VertexTrackProducer.clone(
17+
vertexCompositeCandidates = "KShortEventSelector"
18+
)
1019
LambdaTracks = v0VertexTrackProducer.clone(
11-
vertexCompositeCandidates = "generalV0Candidates:Lambda"
20+
vertexCompositeCandidates = "LambdaEventSelector"
1221
)

0 commit comments

Comments
 (0)