File indexing completed on 2024-09-11 04:33:33
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include "Validation/RecoVertex/interface/V0Validator.h"
0020 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0021
0022 typedef std::vector<TrackingVertex> TrackingVertexCollection;
0023 typedef edm::Ref<TrackingVertexCollection> TrackingVertexRef;
0024 typedef edm::RefVector<edm::HepMCProduct, HepMC::GenVertex> GenVertexRefVector;
0025 typedef edm::RefVector<edm::HepMCProduct, HepMC::GenParticle> GenParticleRefVector;
0026
0027 V0Validator::V0Validator(const edm::ParameterSet& iConfig)
0028 : theDQMRootFileName(iConfig.getUntrackedParameter<std::string>("DQMRootFileName")),
0029 dirName(iConfig.getUntrackedParameter<std::string>("dirName")),
0030 recoRecoToSimCollectionToken_(
0031 consumes<reco::RecoToSimCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
0032 recoSimToRecoCollectionToken_(
0033 consumes<reco::SimToRecoCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
0034 trackingVertexCollection_Token_(
0035 consumes<TrackingVertexCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackingVertexCollection"))),
0036 vec_recoVertex_Token_(
0037 consumes<std::vector<reco::Vertex> >(iConfig.getUntrackedParameter<edm::InputTag>("vertexCollection"))),
0038 recoVertexCompositeCandidateCollection_k0s_Token_(consumes<reco::VertexCompositeCandidateCollection>(
0039 iConfig.getUntrackedParameter<edm::InputTag>("kShortCollection"))),
0040 recoVertexCompositeCandidateCollection_lambda_Token_(consumes<reco::VertexCompositeCandidateCollection>(
0041 iConfig.getUntrackedParameter<edm::InputTag>("lambdaCollection"))) {}
0042
0043 V0Validator::~V0Validator() {}
0044
0045 void V0Validator::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) {
0046 double minKsMass = 0.49767 - 0.07;
0047 double maxKsMass = 0.49767 + 0.07;
0048 double minLamMass = 1.1156 - 0.05;
0049 double maxLamMass = 1.1156 + 0.05;
0050 int ksMassNbins = 100;
0051 double ksMassXmin = minKsMass;
0052 double ksMassXmax = maxKsMass;
0053 int lamMassNbins = 100;
0054 double lamMassXmin = minLamMass;
0055 double lamMassXmax = maxLamMass;
0056
0057 ibooker.cd();
0058 std::string subDirName = V0Validator::dirName + "/K0";
0059 ibooker.setCurrentFolder(subDirName);
0060
0061 candidateEffVsR_num_[V0Validator::KSHORT] =
0062 ibooker.book1D("K0sEffVsR_num", "K^{0}_{S} Efficiency vs #rho", 80, 0., 40.);
0063 candidateEffVsEta_num_[V0Validator::KSHORT] =
0064 ibooker.book1D("K0sEffVsEta_num", "K^{0}_{S} Efficiency vs #eta", 40, -2.5, 2.5);
0065 candidateEffVsPt_num_[V0Validator::KSHORT] =
0066 ibooker.book1D("K0sEffVsPt_num", "K^{0}_{S} Efficiency vs p_{T}", 70, 0., 20.);
0067
0068 candidateTkEffVsR_num_[V0Validator::KSHORT] =
0069 ibooker.book1D("K0sTkEffVsR_num", "K^{0}_{S} Tracking Efficiency vs #rho", 80, 0., 40.);
0070 candidateTkEffVsEta_num_[V0Validator::KSHORT] =
0071 ibooker.book1D("K0sTkEffVsEta_num", "K^{0}_{S} Tracking Efficiency vs #eta", 40, -2.5, 2.5);
0072 candidateTkEffVsPt_num_[V0Validator::KSHORT] =
0073 ibooker.book1D("K0sTkEffVsPt_num", "K^{0}_{S} Tracking Efficiency vs p_{T}", 70, 0., 20.);
0074
0075 candidateEffVsR_denom_[V0Validator::KSHORT] =
0076 ibooker.book1D("K0sEffVsR_denom", "K^{0}_{S} Efficiency vs #rho", 80, 0., 40.);
0077 candidateEffVsEta_denom_[V0Validator::KSHORT] =
0078 ibooker.book1D("K0sEffVsEta_denom", "K^{0}_{S} Efficiency vs #eta", 40, -2.5, 2.5);
0079 candidateEffVsPt_denom_[V0Validator::KSHORT] =
0080 ibooker.book1D("K0sEffVsPt_denom", "K^{0}_{S} Efficiency vs p_{T}", 70, 0., 20.);
0081
0082 candidateFakeVsR_num_[V0Validator::KSHORT] =
0083 ibooker.book1D("K0sFakeVsR_num", "K^{0}_{S} Fake Rate vs #rho", 80, 0., 40.);
0084 candidateFakeVsEta_num_[V0Validator::KSHORT] =
0085 ibooker.book1D("K0sFakeVsEta_num", "K^{0}_{S} Fake Rate vs #eta", 40, -2.5, 2.5);
0086 candidateFakeVsPt_num_[V0Validator::KSHORT] =
0087 ibooker.book1D("K0sFakeVsPt_num", "K^{0}_{S} Fake Rate vs p_{T}", 70, 0., 20.);
0088 candidateTkFakeVsR_num_[V0Validator::KSHORT] =
0089 ibooker.book1D("K0sTkFakeVsR_num", "K^{0}_{S} Tracking Fake Rate vs #rho", 80, 0., 80.);
0090 candidateTkFakeVsEta_num_[V0Validator::KSHORT] =
0091 ibooker.book1D("K0sTkFakeVsEta_num", "K^{0}_{S} Tracking Fake Rate vs #eta", 40, -2.5, 2.5);
0092 candidateTkFakeVsPt_num_[V0Validator::KSHORT] =
0093 ibooker.book1D("K0sTkFakeVsPt_num", "K^{0}_{S} Tracking Fake Rate vs p_{T}", 70, 0., 20.);
0094
0095 candidateFakeVsR_denom_[V0Validator::KSHORT] =
0096 ibooker.book1D("K0sFakeVsR_denom", "K^{0}_{S} Fake Rate vs #rho", 80, 0., 40.);
0097 candidateFakeVsEta_denom_[V0Validator::KSHORT] =
0098 ibooker.book1D("K0sFakeVsEta_denom", "K^{0}_{S} Fake Rate vs #eta", 40, -2.5, 2.5);
0099 candidateFakeVsPt_denom_[V0Validator::KSHORT] =
0100 ibooker.book1D("K0sFakeVsPt_denom", "K^{0}_{S} Fake Rate vs p_{T}", 70, 0., 20.);
0101 nCandidates_[V0Validator::KSHORT] = ibooker.book1D("nK0s", "Number of K^{0}_{S} found per event", 60, 0., 60.);
0102 fakeCandidateMass_[V0Validator::KSHORT] =
0103 ibooker.book1D("ksMassFake", "Mass of fake K0S", ksMassNbins, minKsMass, maxKsMass);
0104 goodCandidateMass[V0Validator::KSHORT] =
0105 ibooker.book1D("ksMassGood", "Mass of good reco K0S", ksMassNbins, minKsMass, maxKsMass);
0106 candidateMassAll[V0Validator::KSHORT] =
0107 ibooker.book1D("ksMassAll", "Invariant mass of all K0S", ksMassNbins, ksMassXmin, ksMassXmax);
0108 candidateFakeDauRadDist_[V0Validator::KSHORT] =
0109 ibooker.book1D("radDistFakeKs", "Production radius of daughter particle of Ks fake", 100, 0., 15.);
0110 candidateStatus_[V0Validator::KSHORT] = ibooker.book1D("ksCandStatus", "Fake type by cand status", 10, 0., 10.);
0111
0112
0113
0114 subDirName = V0Validator::dirName + "/Lambda";
0115 ibooker.setCurrentFolder(subDirName);
0116
0117 candidateEffVsR_num_[V0Validator::LAMBDA] =
0118 ibooker.book1D("LamEffVsR_num", "#Lambda^{0} Efficiency vs #rho", 80, 0., 40.);
0119 candidateEffVsEta_num_[V0Validator::LAMBDA] =
0120 ibooker.book1D("LamEffVsEta_num", "#Lambda^{0} Efficiency vs #eta", 40, -2.5, 2.5);
0121 candidateEffVsPt_num_[V0Validator::LAMBDA] =
0122 ibooker.book1D("LamEffVsPt_num", "#Lambda^{0} Efficiency vs p_{T}", 70, 0., 20.);
0123
0124 candidateTkEffVsR_num_[V0Validator::LAMBDA] =
0125 ibooker.book1D("LamTkEffVsR_num", "#Lambda^{0} TrackingEfficiency vs #rho", 80, 0., 40.);
0126 candidateTkEffVsEta_num_[V0Validator::LAMBDA] =
0127 ibooker.book1D("LamTkEffVsEta_num", "#Lambda^{0} Tracking Efficiency vs #eta", 40, -2.5, 2.5);
0128 candidateTkEffVsPt_num_[V0Validator::LAMBDA] =
0129 ibooker.book1D("LamTkEffVsPt_num", "#Lambda^{0} Tracking Efficiency vs p_{T}", 70, 0., 20.);
0130
0131 candidateEffVsR_denom_[V0Validator::LAMBDA] =
0132 ibooker.book1D("LamEffVsR_denom", "#Lambda^{0} Efficiency vs #rho", 80, 0., 40.);
0133 candidateEffVsEta_denom_[V0Validator::LAMBDA] =
0134 ibooker.book1D("LamEffVsEta_denom", "#Lambda^{0} Efficiency vs #eta", 40, -2.5, 2.5);
0135 candidateEffVsPt_denom_[V0Validator::LAMBDA] =
0136 ibooker.book1D("LamEffVsPt_denom", "#Lambda^{0} Efficiency vs p_{T}", 70, 0., 20.);
0137
0138 candidateFakeVsR_num_[V0Validator::LAMBDA] =
0139 ibooker.book1D("LamFakeVsR_num", "#Lambda^{0} Fake Rate vs #rho", 80, 0., 40.);
0140 candidateFakeVsEta_num_[V0Validator::LAMBDA] =
0141 ibooker.book1D("LamFakeVsEta_num", "#Lambda^{0} Fake Rate vs #eta", 40, -2.5, 2.5);
0142 candidateFakeVsPt_num_[V0Validator::LAMBDA] =
0143 ibooker.book1D("LamFakeVsPt_num", "#Lambda^{0} Fake Rate vs p_{T}", 70, 0., 20.);
0144 candidateTkFakeVsR_num_[V0Validator::LAMBDA] =
0145 ibooker.book1D("LamTkFakeVsR_num", "#Lambda^{0} Tracking Fake Rate vs #rho", 80, 0., 40.);
0146 candidateTkFakeVsEta_num_[V0Validator::LAMBDA] =
0147 ibooker.book1D("LamTkFakeVsEta_num", "#Lambda^{0} Tracking Fake Rate vs #eta", 40, -2.5, 2.5);
0148 candidateTkFakeVsPt_num_[V0Validator::LAMBDA] =
0149 ibooker.book1D("LamTkFakeVsPt_num", "#Lambda^{0} Tracking Fake Rate vs p_{T}", 70, 0., 20.);
0150
0151 candidateFakeVsR_denom_[V0Validator::LAMBDA] =
0152 ibooker.book1D("LamFakeVsR_denom", "#Lambda^{0} Fake Rate vs #rho", 80, 0., 40.);
0153 candidateFakeVsEta_denom_[V0Validator::LAMBDA] =
0154 ibooker.book1D("LamFakeVsEta_denom", "#Lambda^{0} Fake Rate vs #eta", 40, -2.5, 2.5);
0155 candidateFakeVsPt_denom_[V0Validator::LAMBDA] =
0156 ibooker.book1D("LamFakeVsPt_denom", "#Lambda^{0} Fake Rate vs p_{T}", 70, 0., 20.);
0157
0158 nCandidates_[V0Validator::LAMBDA] = ibooker.book1D("nLam", "Number of #Lambda^{0} found per event", 60, 0., 60.);
0159 fakeCandidateMass_[V0Validator::LAMBDA] =
0160 ibooker.book1D("lamMassFake", "Mass of fake Lambda", lamMassNbins, minLamMass, maxLamMass);
0161 goodCandidateMass[V0Validator::LAMBDA] =
0162 ibooker.book1D("lamMassGood", "Mass of good Lambda", lamMassNbins, minLamMass, maxLamMass);
0163
0164 candidateMassAll[V0Validator::LAMBDA] =
0165 ibooker.book1D("lamMassAll", "Invariant mass of all #Lambda^{0}", lamMassNbins, lamMassXmin, lamMassXmax);
0166 candidateFakeDauRadDist_[V0Validator::LAMBDA] =
0167 ibooker.book1D("radDistFakeLam", "Production radius of daughter particle of Lam fake", 100, 0., 15.);
0168
0169 candidateStatus_[V0Validator::LAMBDA] = ibooker.book1D("ksCandStatus", "Fake type by cand status", 10, 0., 10.);
0170 }
0171
0172 void V0Validator::doFakeRates(const reco::VertexCompositeCandidateCollection& collection,
0173 const reco::RecoToSimCollection& recotosimCollection,
0174 V0Type v0_type,
0175 int particle_pdgid,
0176 int misreconstructed_particle_pdgid) {
0177 using namespace edm;
0178
0179 int numCandidateFound = 0;
0180 double mass = 0.;
0181 float CandidatepT = 0.;
0182 float CandidateEta = 0.;
0183 float CandidateR = 0.;
0184 int CandidateStatus = 0;
0185 const unsigned int NUM_DAUGHTERS = 2;
0186 if (!collection.empty()) {
0187 for (reco::VertexCompositeCandidateCollection::const_iterator iCandidate = collection.begin();
0188 iCandidate != collection.end();
0189 iCandidate++) {
0190
0191 mass = iCandidate->mass();
0192 CandidatepT = (sqrt(iCandidate->momentum().perp2()));
0193 CandidateEta = iCandidate->momentum().eta();
0194 CandidateR = (sqrt(iCandidate->vertex().perp2()));
0195 candidateMassAll[v0_type]->Fill(mass);
0196 CandidateStatus = 0;
0197
0198 std::array<reco::TrackRef, NUM_DAUGHTERS> theDaughterTracks = {
0199 {(*(dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(0)))).track(),
0200 (*(dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(1)))).track()}};
0201
0202 TrackingParticleRef tpref;
0203 TrackingParticleRef firstDauTP;
0204 TrackingVertexRef candidateVtx;
0205
0206 std::array<double, NUM_DAUGHTERS> radDist;
0207
0208 for (View<reco::Track>::size_type i = 0; i < theDaughterTracks.size(); ++i) {
0209 radDist = {{-1., -1.}};
0210
0211 RefToBase<reco::Track> track(theDaughterTracks.at(i));
0212
0213 if (recotosimCollection.find(track) != recotosimCollection.end()) {
0214 const std::vector<std::pair<TrackingParticleRef, double> >& tp = recotosimCollection[track];
0215 if (!tp.empty()) {
0216 tpref = tp.begin()->first;
0217
0218 TrackingVertexRef parentVertex = tpref->parentVertex();
0219 if (parentVertex.isNonnull()) {
0220 radDist[i] = parentVertex->position().R();
0221 if (candidateVtx.isNonnull()) {
0222 if (candidateVtx->position() == parentVertex->position()) {
0223 if (parentVertex->nDaughterTracks() == 2) {
0224 if (parentVertex->nSourceTracks() == 0) {
0225
0226
0227
0228 CandidateStatus = 6;
0229 }
0230
0231 for (TrackingVertex::tp_iterator iTP = parentVertex->sourceTracks_begin();
0232 iTP != parentVertex->sourceTracks_end();
0233 iTP++) {
0234 if (abs((*iTP)->pdgId()) == particle_pdgid) {
0235 CandidateStatus = 1;
0236 numCandidateFound += 1.;
0237 goodCandidateMass[v0_type]->Fill(mass);
0238 } else {
0239 CandidateStatus = 2;
0240 if (abs((*iTP)->pdgId()) == misreconstructed_particle_pdgid) {
0241 CandidateStatus = 7;
0242 }
0243 }
0244 }
0245 } else {
0246
0247
0248 CandidateStatus = 3;
0249 }
0250 } else {
0251
0252
0253 CandidateStatus = 4;
0254 }
0255 } else {
0256
0257
0258
0259 candidateVtx = parentVertex;
0260 firstDauTP = tpref;
0261 }
0262 }
0263 }
0264 } else {
0265 CandidateStatus = 5;
0266 }
0267 }
0268
0269
0270 if (CandidateStatus > 1) {
0271 candidateFakeVsR_num_[v0_type]->Fill(CandidateR);
0272 candidateFakeVsEta_num_[v0_type]->Fill(CandidateEta);
0273 candidateFakeVsPt_num_[v0_type]->Fill(CandidatepT);
0274 candidateStatus_[v0_type]->Fill((float)CandidateStatus);
0275 fakeCandidateMass_[v0_type]->Fill(mass);
0276 for (auto distance : radDist) {
0277 if (distance > 0)
0278 candidateFakeDauRadDist_[v0_type]->Fill(distance);
0279 }
0280 }
0281 if (CandidateStatus == 5) {
0282 candidateTkFakeVsR_num_[v0_type]->Fill(CandidateR);
0283 candidateTkFakeVsEta_num_[v0_type]->Fill(CandidateEta);
0284 candidateTkFakeVsPt_num_[v0_type]->Fill(CandidatepT);
0285 }
0286 candidateFakeVsR_denom_[v0_type]->Fill(CandidateR);
0287 candidateFakeVsEta_denom_[v0_type]->Fill(CandidateEta);
0288 candidateFakeVsPt_denom_[v0_type]->Fill(CandidatepT);
0289 }
0290 }
0291 nCandidates_[v0_type]->Fill((float)numCandidateFound);
0292 }
0293
0294 void V0Validator::doEfficiencies(const TrackingVertexCollection& gen_vertices,
0295 V0Type v0_type,
0296 int parent_particle_id,
0297 int first_daughter_id,
0298 int second_daughter_id,
0299 const reco::VertexCompositeCandidateCollection& collection,
0300 const reco::SimToRecoCollection& simtorecoCollection) {
0301
0302
0303
0304
0305
0306
0307 std::set<V0Couple> reconstructed_V0_couples;
0308 if (!collection.empty()) {
0309 for (reco::VertexCompositeCandidateCollection::const_iterator iCandidate = collection.begin();
0310 iCandidate != collection.end();
0311 iCandidate++) {
0312 reconstructed_V0_couples.insert(
0313 V0Couple((dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(0)))->track(),
0314 (dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(1)))->track()));
0315 }
0316 }
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328 unsigned int candidateEff[2] = {0, 0};
0329 for (auto const& gen_vertex : gen_vertices) {
0330 if (gen_vertex.eventId().bunchCrossing() != 0)
0331 continue;
0332 if (gen_vertex.nDaughterTracks() != 2)
0333 continue;
0334 for (TrackingVertex::tp_iterator source = gen_vertex.sourceTracks_begin(); source != gen_vertex.sourceTracks_end();
0335 ++source) {
0336 if (std::abs((*source)->pdgId()) == parent_particle_id) {
0337 if ((std::abs((gen_vertex.daughterTracks().at(0))->pdgId()) == first_daughter_id &&
0338 std::abs((gen_vertex.daughterTracks().at(1))->pdgId()) == second_daughter_id) ||
0339 (std::abs((gen_vertex.daughterTracks().at(0))->pdgId()) == second_daughter_id &&
0340 std::abs((gen_vertex.daughterTracks().at(1))->pdgId()) == first_daughter_id)) {
0341 if ((std::abs((gen_vertex.daughterTracks().at(0))->momentum().eta()) < 2.4 &&
0342 gen_vertex.daughterTracks().at(0)->pt() > 0.9) &&
0343 (std::abs((gen_vertex.daughterTracks().at(1))->momentum().eta()) < 2.4 &&
0344 gen_vertex.daughterTracks().at(1)->pt() > 0.9)) {
0345
0346 float candidateGenpT = sqrt((*source)->momentum().perp2());
0347 float candidateGenEta = (*source)->momentum().eta();
0348 float candidateGenR = sqrt((*source)->vertex().perp2());
0349 candidateEffVsPt_denom_[v0_type]->Fill(candidateGenpT);
0350 candidateEffVsEta_denom_[v0_type]->Fill(candidateGenEta);
0351 candidateEffVsR_denom_[v0_type]->Fill(candidateGenR);
0352
0353 std::array<reco::TrackRef, 2> reco_daughter;
0354
0355 for (unsigned int daughter = 0; daughter < 2; ++daughter) {
0356 if (simtorecoCollection.find(gen_vertex.daughterTracks()[daughter]) != simtorecoCollection.end()) {
0357 if (!simtorecoCollection[gen_vertex.daughterTracks()[daughter]].empty()) {
0358 candidateEff[daughter] = 1;
0359 reco_daughter[daughter] =
0360 simtorecoCollection[gen_vertex.daughterTracks()[daughter]].begin()->first.castTo<reco::TrackRef>();
0361 }
0362 } else {
0363 candidateEff[daughter] = 2;
0364 }
0365 }
0366 if ((candidateEff[0] == 1 && candidateEff[1] == 1) && (reco_daughter[0].key() != reco_daughter[1].key()) &&
0367 (reconstructed_V0_couples.find(V0Couple(reco_daughter[0], reco_daughter[1])) !=
0368 reconstructed_V0_couples.end())) {
0369 candidateEffVsPt_num_[v0_type]->Fill(candidateGenpT);
0370 candidateEffVsEta_num_[v0_type]->Fill(candidateGenEta);
0371 candidateEffVsR_num_[v0_type]->Fill(candidateGenR);
0372 }
0373 }
0374 }
0375 }
0376 }
0377 }
0378 }
0379
0380 void V0Validator::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0381 using std::cout;
0382 using std::endl;
0383 using namespace edm;
0384 using namespace std;
0385
0386
0387 Handle<reco::RecoToSimCollection> recotosimCollectionH;
0388 iEvent.getByToken(recoRecoToSimCollectionToken_, recotosimCollectionH);
0389
0390 Handle<reco::SimToRecoCollection> simtorecoCollectionH;
0391 iEvent.getByToken(recoSimToRecoCollectionToken_, simtorecoCollectionH);
0392
0393
0394 edm::Handle<TrackingVertexCollection> TVCollectionH;
0395 iEvent.getByToken(trackingVertexCollection_Token_, TVCollectionH);
0396
0397
0398 edm::Handle<std::vector<reco::Vertex> > primaryVtxCollectionH;
0399 iEvent.getByToken(vec_recoVertex_Token_, primaryVtxCollectionH);
0400
0401 std::vector<reco::Vertex>::const_iterator iVtxPH = primaryVtxCollectionH->begin();
0402 for (std::vector<reco::Vertex>::const_iterator iVtx = primaryVtxCollectionH->begin();
0403 iVtx < primaryVtxCollectionH->end();
0404 iVtx++) {
0405 if (primaryVtxCollectionH->size() > 1) {
0406 if (iVtx->tracksSize() > iVtxPH->tracksSize()) {
0407 iVtxPH = iVtx;
0408 }
0409 } else
0410 iVtxPH = iVtx;
0411 }
0412
0413
0414 edm::Handle<reco::VertexCompositeCandidateCollection> k0sCollection;
0415 edm::Handle<reco::VertexCompositeCandidateCollection> lambdaCollection;
0416 iEvent.getByToken(recoVertexCompositeCandidateCollection_k0s_Token_, k0sCollection);
0417 iEvent.getByToken(recoVertexCompositeCandidateCollection_lambda_Token_, lambdaCollection);
0418
0419
0420
0421
0422
0423 if (k0sCollection.isValid()) {
0424 doFakeRates(*k0sCollection.product(), *recotosimCollectionH.product(), V0Type::KSHORT, 310, 3122);
0425 doEfficiencies(*TVCollectionH.product(),
0426 V0Type::KSHORT,
0427 310,
0428 211,
0429 211,
0430 *k0sCollection.product(),
0431 *simtorecoCollectionH.product());
0432 }
0433 if (lambdaCollection.isValid()) {
0434 doFakeRates(*lambdaCollection.product(), *recotosimCollectionH.product(), V0Type::LAMBDA, 3122, 310);
0435 doEfficiencies(*TVCollectionH.product(),
0436 V0Type::LAMBDA,
0437 3122,
0438 211,
0439 2212,
0440 *lambdaCollection.product(),
0441 *simtorecoCollectionH.product());
0442 }
0443 }
0444
0445
0446