Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:10

0001 // -*- C++ -*-
0002 //
0003 // Package:    DataFormats/Scouting
0004 // Class:      TestReadRun2Scouting
0005 //
0006 /**\class edmtest::TestReadRun2Scouting
0007   Description: Used as part of tests that ensure the run 2 Scouting
0008   data formats can be persistently written and in a subsequent process
0009   read. First, this is done using the current release version for writing
0010   and reading. In addition, the output file of the write process should
0011   be saved permanently each time a run 2 Scouting persistent data
0012   format changes. In unit tests, we read each of those saved files to verify
0013   that the current releases can read older versions of these data formats.
0014 */
0015 // Original Author:  W. David Dagenhart
0016 //         Created:  2 June 2023
0017 
0018 #include "DataFormats/Scouting/interface/ScoutingCaloJet.h"
0019 #include "DataFormats/Scouting/interface/ScoutingElectron.h"
0020 #include "DataFormats/Scouting/interface/ScoutingMuon.h"
0021 #include "DataFormats/Scouting/interface/ScoutingParticle.h"
0022 #include "DataFormats/Scouting/interface/ScoutingPFJet.h"
0023 #include "DataFormats/Scouting/interface/ScoutingPhoton.h"
0024 #include "DataFormats/Scouting/interface/ScoutingTrack.h"
0025 #include "DataFormats/Scouting/interface/ScoutingVertex.h"
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0033 #include "FWCore/Utilities/interface/EDGetToken.h"
0034 #include "FWCore/Utilities/interface/Exception.h"
0035 #include "FWCore/Utilities/interface/InputTag.h"
0036 #include "FWCore/Utilities/interface/StreamID.h"
0037 
0038 #include <vector>
0039 
0040 namespace edmtest {
0041 
0042   class TestReadRun2Scouting : public edm::global::EDAnalyzer<> {
0043   public:
0044     TestReadRun2Scouting(edm::ParameterSet const&);
0045     void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;
0046     static void fillDescriptions(edm::ConfigurationDescriptions&);
0047 
0048   private:
0049     void analyzeCaloJets(edm::Event const&) const;
0050     void analyzeElectrons(edm::Event const&) const;
0051     void analyzeMuons(edm::Event const&) const;
0052     void analyzeParticles(edm::Event const&) const;
0053     void analyzePFJets(edm::Event const&) const;
0054     void analyzePhotons(edm::Event const&) const;
0055     void analyzeTracks(edm::Event const&) const;
0056     void analyzeVertexes(edm::Event const&) const;
0057 
0058     void throwWithMessageFromConstructor(const char*) const;
0059     void throwWithMessage(const char*) const;
0060 
0061     // These expected values are meaningless other than we use them
0062     // to check that values read from persistent storage match the values
0063     // we know were written.
0064 
0065     const std::vector<double> expectedCaloJetsValues_;
0066     const edm::EDGetTokenT<std::vector<ScoutingCaloJet>> caloJetsToken_;
0067 
0068     const std::vector<double> expectedElectronFloatingPointValues_;
0069     const std::vector<int> expectedElectronIntegralValues_;
0070     const edm::EDGetTokenT<std::vector<ScoutingElectron>> electronsToken_;
0071 
0072     const int inputMuonClassVersion_;
0073     const std::vector<double> expectedMuonFloatingPointValues_;
0074     const std::vector<int> expectedMuonIntegralValues_;
0075     const edm::EDGetTokenT<std::vector<ScoutingMuon>> muonsToken_;
0076 
0077     const std::vector<double> expectedParticleFloatingPointValues_;
0078     const std::vector<int> expectedParticleIntegralValues_;
0079     const edm::EDGetTokenT<std::vector<ScoutingParticle>> particlesToken_;
0080 
0081     const std::vector<double> expectedPFJetFloatingPointValues_;
0082     const std::vector<int> expectedPFJetIntegralValues_;
0083     const edm::EDGetTokenT<std::vector<ScoutingPFJet>> pfJetsToken_;
0084 
0085     const std::vector<double> expectedPhotonFloatingPointValues_;
0086     const edm::EDGetTokenT<std::vector<ScoutingPhoton>> photonsToken_;
0087 
0088     const int inputTrackClassVersion_;
0089     const std::vector<double> expectedTrackFloatingPointValues_;
0090     const std::vector<int> expectedTrackIntegralValues_;
0091     const edm::EDGetTokenT<std::vector<ScoutingTrack>> tracksToken_;
0092 
0093     const int inputVertexClassVersion_;
0094     const std::vector<double> expectedVertexFloatingPointValues_;
0095     const std::vector<int> expectedVertexIntegralValues_;
0096     const edm::EDGetTokenT<std::vector<ScoutingVertex>> vertexesToken_;
0097   };
0098 
0099   TestReadRun2Scouting::TestReadRun2Scouting(edm::ParameterSet const& iPSet)
0100       : expectedCaloJetsValues_(iPSet.getParameter<std::vector<double>>("expectedCaloJetsValues")),
0101         caloJetsToken_(consumes(iPSet.getParameter<edm::InputTag>("caloJetsTag"))),
0102         expectedElectronFloatingPointValues_(
0103             iPSet.getParameter<std::vector<double>>("expectedElectronFloatingPointValues")),
0104         expectedElectronIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedElectronIntegralValues")),
0105         electronsToken_(consumes(iPSet.getParameter<edm::InputTag>("electronsTag"))),
0106         inputMuonClassVersion_(iPSet.getParameter<int>("muonClassVersion")),
0107         expectedMuonFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedMuonFloatingPointValues")),
0108         expectedMuonIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedMuonIntegralValues")),
0109         muonsToken_(consumes(iPSet.getParameter<edm::InputTag>("muonsTag"))),
0110         expectedParticleFloatingPointValues_(
0111             iPSet.getParameter<std::vector<double>>("expectedParticleFloatingPointValues")),
0112         expectedParticleIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedParticleIntegralValues")),
0113         particlesToken_(consumes(iPSet.getParameter<edm::InputTag>("particlesTag"))),
0114         expectedPFJetFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedPFJetFloatingPointValues")),
0115         expectedPFJetIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedPFJetIntegralValues")),
0116         pfJetsToken_(consumes(iPSet.getParameter<edm::InputTag>("pfJetsTag"))),
0117         expectedPhotonFloatingPointValues_(
0118             iPSet.getParameter<std::vector<double>>("expectedPhotonFloatingPointValues")),
0119         photonsToken_(consumes(iPSet.getParameter<edm::InputTag>("photonsTag"))),
0120         inputTrackClassVersion_(iPSet.getParameter<int>("trackClassVersion")),
0121         expectedTrackFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedTrackFloatingPointValues")),
0122         expectedTrackIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedTrackIntegralValues")),
0123         tracksToken_(consumes(iPSet.getParameter<edm::InputTag>("tracksTag"))),
0124         inputVertexClassVersion_(iPSet.getParameter<int>("vertexClassVersion")),
0125         expectedVertexFloatingPointValues_(
0126             iPSet.getParameter<std::vector<double>>("expectedVertexFloatingPointValues")),
0127         expectedVertexIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedVertexIntegralValues")),
0128         vertexesToken_(consumes(iPSet.getParameter<edm::InputTag>("vertexesTag"))) {
0129     if (expectedCaloJetsValues_.size() != 16) {
0130       throwWithMessageFromConstructor("test configuration error, expectedCaloJetsValues must have size 16");
0131     }
0132     if (expectedElectronFloatingPointValues_.size() != 14) {
0133       throwWithMessageFromConstructor(
0134           "test configuration error, expectedElectronFloatingPointValues must have size 14");
0135     }
0136     if (expectedElectronIntegralValues_.size() != 2) {
0137       throwWithMessageFromConstructor("test configuration error, expectedElectronIntegralValues must have size 2");
0138     }
0139     if (expectedMuonFloatingPointValues_.size() != 23) {
0140       throwWithMessageFromConstructor("test configuration error, expectedMuonFloatingPointValues must have size 23");
0141     }
0142     if (expectedMuonIntegralValues_.size() != 8) {
0143       throwWithMessageFromConstructor("test configuration error, expectedMuonIntegralValues must have size 8");
0144     }
0145     if (expectedParticleFloatingPointValues_.size() != 4) {
0146       throwWithMessageFromConstructor("test configuration error, expectedParticleFloatingPointValues must have size 4");
0147     }
0148     if (expectedParticleIntegralValues_.size() != 2) {
0149       throwWithMessageFromConstructor("test configuration error, expectedParticleIntegralValues must have size 2");
0150     }
0151     if (expectedPFJetFloatingPointValues_.size() != 15) {
0152       throwWithMessageFromConstructor("test configuration error, expectedPFJetFloatingPointValues must have size 15");
0153     }
0154     if (expectedPFJetIntegralValues_.size() != 8) {
0155       throwWithMessageFromConstructor("test configuration error, expectedPFJetIntegralValues must have size 8");
0156     }
0157     if (expectedPhotonFloatingPointValues_.size() != 8) {
0158       throwWithMessageFromConstructor("test configuration error, expectedPhotonFloatingPointValues must have size 8");
0159     }
0160     if (expectedTrackFloatingPointValues_.size() != 16) {
0161       throwWithMessageFromConstructor("test configuration error, expectedTrackFloatingPointValues must have size 16");
0162     }
0163     if (expectedTrackIntegralValues_.size() != 4) {
0164       throwWithMessageFromConstructor("test configuration error, expectedTrackIntegralValues must have size 4");
0165     }
0166     if (expectedVertexFloatingPointValues_.size() != 7) {
0167       throwWithMessageFromConstructor("test configuration error, expectedVertexFloatingPointValues must have size 7");
0168     }
0169     if (expectedVertexIntegralValues_.size() != 3) {
0170       throwWithMessageFromConstructor("test configuration error, expectedVertexIntegralValues must have size 3");
0171     }
0172   }
0173 
0174   void TestReadRun2Scouting::analyze(edm::StreamID, edm::Event const& iEvent, edm::EventSetup const&) const {
0175     analyzeCaloJets(iEvent);
0176     analyzeElectrons(iEvent);
0177     analyzeMuons(iEvent);
0178     analyzeParticles(iEvent);
0179     analyzePFJets(iEvent);
0180     analyzePhotons(iEvent);
0181     analyzeTracks(iEvent);
0182     analyzeVertexes(iEvent);
0183   }
0184 
0185   void TestReadRun2Scouting::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0186     edm::ParameterSetDescription desc;
0187     desc.add<std::vector<double>>("expectedCaloJetsValues");
0188     desc.add<edm::InputTag>("caloJetsTag");
0189     desc.add<std::vector<double>>("expectedElectronFloatingPointValues");
0190     desc.add<std::vector<int>>("expectedElectronIntegralValues");
0191     desc.add<edm::InputTag>("electronsTag");
0192     desc.add<int>("muonClassVersion");
0193     desc.add<std::vector<double>>("expectedMuonFloatingPointValues");
0194     desc.add<std::vector<int>>("expectedMuonIntegralValues");
0195     desc.add<edm::InputTag>("muonsTag");
0196     desc.add<std::vector<double>>("expectedParticleFloatingPointValues");
0197     desc.add<std::vector<int>>("expectedParticleIntegralValues");
0198     desc.add<edm::InputTag>("particlesTag");
0199     desc.add<std::vector<double>>("expectedPFJetFloatingPointValues");
0200     desc.add<std::vector<int>>("expectedPFJetIntegralValues");
0201     desc.add<edm::InputTag>("pfJetsTag");
0202     desc.add<std::vector<double>>("expectedPhotonFloatingPointValues");
0203     desc.add<edm::InputTag>("photonsTag");
0204     desc.add<int>("trackClassVersion");
0205     desc.add<std::vector<double>>("expectedTrackFloatingPointValues");
0206     desc.add<std::vector<int>>("expectedTrackIntegralValues");
0207     desc.add<edm::InputTag>("tracksTag");
0208     desc.add<int>("vertexClassVersion");
0209     desc.add<std::vector<double>>("expectedVertexFloatingPointValues");
0210     desc.add<std::vector<int>>("expectedVertexIntegralValues");
0211     desc.add<edm::InputTag>("vertexesTag");
0212     descriptions.addDefault(desc);
0213   }
0214 
0215   void TestReadRun2Scouting::analyzeCaloJets(edm::Event const& iEvent) const {
0216     auto const& caloJets = iEvent.get(caloJetsToken_);
0217     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0218     if (caloJets.size() != vectorSize) {
0219       throwWithMessage("analyzeCaloJets, caloJets does not have expected size");
0220     }
0221     unsigned int i = 0;
0222     for (auto const& caloJet : caloJets) {
0223       double offset = static_cast<double>(iEvent.id().event() + i);
0224 
0225       if (caloJet.pt() != expectedCaloJetsValues_[0] + offset) {
0226         throwWithMessage("analyzeCaloJets, pt does not equal expected value");
0227       }
0228       if (caloJet.eta() != expectedCaloJetsValues_[1] + offset) {
0229         throwWithMessage("analyzeCaloJets, eta does not equal expected value");
0230       }
0231       if (caloJet.phi() != expectedCaloJetsValues_[2] + offset) {
0232         throwWithMessage("analyzeCaloJets, phi does not equal expected value");
0233       }
0234       if (caloJet.m() != expectedCaloJetsValues_[3] + offset) {
0235         throwWithMessage("analyzeCaloJets, m does not equal expected value");
0236       }
0237       if (caloJet.jetArea() != expectedCaloJetsValues_[4] + offset) {
0238         throwWithMessage("analyzeCaloJets, jetArea does not equal expected value");
0239       }
0240       if (caloJet.maxEInEmTowers() != expectedCaloJetsValues_[5] + offset) {
0241         throwWithMessage("analyzeCaloJets,  maxEInEmTowers() does not equal expected value");
0242       }
0243       if (caloJet.maxEInHadTowers() != expectedCaloJetsValues_[6] + offset) {
0244         throwWithMessage("analyzeCaloJets,  maxEInHadTowers does not equal expected value");
0245       }
0246       if (caloJet.hadEnergyInHB() != expectedCaloJetsValues_[7] + offset) {
0247         throwWithMessage("analyzeCaloJets, hadEnergyInHB does not equal expected value");
0248       }
0249       if (caloJet.hadEnergyInHE() != expectedCaloJetsValues_[8] + offset) {
0250         throwWithMessage("analyzeCaloJets, hadEnergyInHE does not equal expected value");
0251       }
0252       if (caloJet.hadEnergyInHF() != expectedCaloJetsValues_[9] + offset) {
0253         throwWithMessage("analyzeCaloJets, hadEnergyInHF does not equal expected value");
0254       }
0255       if (caloJet.emEnergyInEB() != expectedCaloJetsValues_[10] + offset) {
0256         throwWithMessage("analyzeCaloJets, emEnergyInEB does not equal expected value");
0257       }
0258       if (caloJet.emEnergyInEE() != expectedCaloJetsValues_[11] + offset) {
0259         throwWithMessage("analyzeCaloJets, emEnergyInEE does not equal expected value");
0260       }
0261       if (caloJet.emEnergyInHF() != expectedCaloJetsValues_[12] + offset) {
0262         throwWithMessage("analyzeCaloJets, emEnergyInHF does not equal expected value");
0263       }
0264       if (caloJet.towersArea() != expectedCaloJetsValues_[13] + offset) {
0265         throwWithMessage("analyzeCaloJets, towersArea does not equal expected value");
0266       }
0267       if (caloJet.mvaDiscriminator() != expectedCaloJetsValues_[14] + offset) {
0268         throwWithMessage("analyzeCaloJets,  mvaDiscriminator does not equal expected value");
0269       }
0270       if (caloJet.btagDiscriminator() != expectedCaloJetsValues_[15] + offset) {
0271         throwWithMessage("analyzeCaloJets,  btagDiscriminator does not equal expected value");
0272       }
0273       ++i;
0274     }
0275   }
0276 
0277   void TestReadRun2Scouting::analyzeElectrons(edm::Event const& iEvent) const {
0278     auto const& electrons = iEvent.get(electronsToken_);
0279     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0280     if (electrons.size() != vectorSize) {
0281       throwWithMessage("analyzeElectrons, electrons does not have expected size");
0282     }
0283     unsigned int i = 0;
0284     for (auto const& electron : electrons) {
0285       double offset = static_cast<double>(iEvent.id().event() + i);
0286       int iOffset = static_cast<int>(iEvent.id().event() + i);
0287 
0288       if (electron.pt() != expectedElectronFloatingPointValues_[0] + offset) {
0289         throwWithMessage("analyzeElectrons, pt does not equal expected value");
0290       }
0291       if (electron.eta() != expectedElectronFloatingPointValues_[1] + offset) {
0292         throwWithMessage("analyzeElectrons, eta does not equal expected value");
0293       }
0294       if (electron.phi() != expectedElectronFloatingPointValues_[2] + offset) {
0295         throwWithMessage("analyzeElectrons, phi does not equal expected value");
0296       }
0297       if (electron.m() != expectedElectronFloatingPointValues_[3] + offset) {
0298         throwWithMessage("analyzeElectrons, m does not equal expected value");
0299       }
0300       if (electron.d0() != expectedElectronFloatingPointValues_[4] + offset) {
0301         throwWithMessage("analyzeElectrons, d0 does not equal expected value");
0302       }
0303       if (electron.dz() != expectedElectronFloatingPointValues_[5] + offset) {
0304         throwWithMessage("analyzeElectrons, dz does not equal expected value");
0305       }
0306       if (electron.dEtaIn() != expectedElectronFloatingPointValues_[6] + offset) {
0307         throwWithMessage("analyzeElectrons,  dEtaIn does not equal expected value");
0308       }
0309       if (electron.dPhiIn() != expectedElectronFloatingPointValues_[7] + offset) {
0310         throwWithMessage("analyzeElectrons, dPhiIn does not equal expected value");
0311       }
0312       if (electron.sigmaIetaIeta() != expectedElectronFloatingPointValues_[8] + offset) {
0313         throwWithMessage("analyzeElectrons, sigmaIetaIeta does not equal expected value");
0314       }
0315       if (electron.hOverE() != expectedElectronFloatingPointValues_[9] + offset) {
0316         throwWithMessage("analyzeElectrons, hOverE does not equal expected value");
0317       }
0318       if (electron.ooEMOop() != expectedElectronFloatingPointValues_[10] + offset) {
0319         throwWithMessage("analyzeElectrons, ooEMOop does not equal expected value");
0320       }
0321       if (electron.missingHits() != expectedElectronIntegralValues_[0] + iOffset) {
0322         throwWithMessage("analyzeElectrons, missingHits does not equal expected value");
0323       }
0324       if (electron.charge() != expectedElectronIntegralValues_[1] + iOffset) {
0325         throwWithMessage("analyzeElectrons, charge does not equal expected value");
0326       }
0327       if (electron.ecalIso() != expectedElectronFloatingPointValues_[11] + offset) {
0328         throwWithMessage("analyzeElectrons, ecalIso does not equal expected value");
0329       }
0330       if (electron.hcalIso() != expectedElectronFloatingPointValues_[12] + offset) {
0331         throwWithMessage("analyzeElectrons, hcalIso does not equal expected value");
0332       }
0333       if (electron.trackIso() != expectedElectronFloatingPointValues_[13] + offset) {
0334         throwWithMessage("analyzeElectrons, trackIso does not equal expected value");
0335       }
0336       ++i;
0337     }
0338   }
0339 
0340   void TestReadRun2Scouting::analyzeMuons(edm::Event const& iEvent) const {
0341     auto const& muons = iEvent.get(muonsToken_);
0342     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0343     if (muons.size() != vectorSize) {
0344       throwWithMessage("analyzeMuons, muons does not have expected size");
0345     }
0346     unsigned int i = 0;
0347     for (auto const& muon : muons) {
0348       double offset = static_cast<double>(iEvent.id().event() + i);
0349       int iOffset = static_cast<int>(iEvent.id().event() + i);
0350 
0351       if (muon.pt() != expectedMuonFloatingPointValues_[0] + offset) {
0352         throwWithMessage("analyzeMuons, pt does not equal expected value");
0353       }
0354       if (muon.eta() != expectedMuonFloatingPointValues_[1] + offset) {
0355         throwWithMessage("analyzeMuons, eta does not equal expected value");
0356       }
0357       if (muon.phi() != expectedMuonFloatingPointValues_[2] + offset) {
0358         throwWithMessage("analyzeMuons, phi does not equal expected value");
0359       }
0360       if (muon.m() != expectedMuonFloatingPointValues_[3] + offset) {
0361         throwWithMessage("analyzeMuons, m does not equal expected value");
0362       }
0363       if (muon.ecalIso() != expectedMuonFloatingPointValues_[4] + offset) {
0364         throwWithMessage("analyzeMuons, ecalIso does not equal expected value");
0365       }
0366       if (muon.hcalIso() != expectedMuonFloatingPointValues_[5] + offset) {
0367         throwWithMessage("analyzeMuons, hcalIso does not equal expected value");
0368       }
0369       if (muon.trackIso() != expectedMuonFloatingPointValues_[6] + offset) {
0370         throwWithMessage("analyzeMuons, trackIso does not equal expected value");
0371       }
0372       if (muon.chi2() != expectedMuonFloatingPointValues_[7] + offset) {
0373         throwWithMessage("analyzeMuons, chi2 does not equal expected value");
0374       }
0375       if (muon.ndof() != expectedMuonFloatingPointValues_[8] + offset) {
0376         throwWithMessage("analyzeMuons, ndof does not equal expected value");
0377       }
0378       if (muon.charge() != expectedMuonIntegralValues_[0] + iOffset) {
0379         throwWithMessage("analyzeMuons, charge does not equal expected value");
0380       }
0381       if (muon.dxy() != expectedMuonFloatingPointValues_[9] + offset) {
0382         throwWithMessage("analyzeMuons, dxy does not equal expected value");
0383       }
0384       if (muon.dz() != expectedMuonFloatingPointValues_[10] + offset) {
0385         throwWithMessage("analyzeMuons, dz does not equal expected value");
0386       }
0387       if (muon.nValidMuonHits() != expectedMuonIntegralValues_[1] + iOffset) {
0388         throwWithMessage("analyzeMuons, nValidMuonHits does not equal expected value");
0389       }
0390       if (muon.nValidPixelHits() != expectedMuonIntegralValues_[2] + iOffset) {
0391         throwWithMessage("analyzeMuons, nValidPixelHits does not equal expected value");
0392       }
0393       if (muon.nMatchedStations() != expectedMuonIntegralValues_[3] + iOffset) {
0394         throwWithMessage("analyzeMuons, nMatchedStations does not equal expected value");
0395       }
0396       if (muon.nTrackerLayersWithMeasurement() != expectedMuonIntegralValues_[4] + iOffset) {
0397         throwWithMessage("analyzeMuons, nTrackerLayersWithMeasurement does not equal expected value");
0398       }
0399       if (muon.type() != expectedMuonIntegralValues_[5] + iOffset) {
0400         throwWithMessage("analyzeMuons, type does not equal expected value");
0401       }
0402 
0403       if (inputMuonClassVersion_ > 2) {
0404         if (muon.nValidStripHits() != expectedMuonIntegralValues_[6] + iOffset) {
0405           throwWithMessage("analyzeMuons, nValidStripHits does not equal expected value");
0406         }
0407         if (muon.trk_qoverp() != expectedMuonFloatingPointValues_[11] + offset) {
0408           throwWithMessage("analyzeMuons, trk_qoverp does not equal expected value");
0409         }
0410         if (muon.trk_lambda() != expectedMuonFloatingPointValues_[12] + offset) {
0411           throwWithMessage("analyzeMuons, trk_lambda does not equal expected value");
0412         }
0413         if (muon.trk_pt() != expectedMuonFloatingPointValues_[13] + offset) {
0414           throwWithMessage("analyzeMuons, trk_pt does not equal expected value");
0415         }
0416         if (muon.trk_phi() != expectedMuonFloatingPointValues_[14] + offset) {
0417           throwWithMessage("analyzeMuons, trk_phi does not equal expected value");
0418         }
0419         if (muon.trk_eta() != expectedMuonFloatingPointValues_[15] + offset) {
0420           throwWithMessage("analyzeMuons, trk_eta does not equal expected value");
0421         }
0422         if (muon.dxyError() != expectedMuonFloatingPointValues_[16] + offset) {
0423           throwWithMessage("analyzeMuons, dxyError does not equal expected value");
0424         }
0425         if (muon.dzError() != expectedMuonFloatingPointValues_[17] + offset) {
0426           throwWithMessage("analyzeMuons, dzError does not equal expected value");
0427         }
0428         if (muon.trk_qoverpError() != expectedMuonFloatingPointValues_[18] + offset) {
0429           throwWithMessage("analyzeMuons, trk_qoverpError does not equal expected value");
0430         }
0431         if (muon.trk_lambdaError() != expectedMuonFloatingPointValues_[19] + offset) {
0432           throwWithMessage("analyzeMuons, trk_lambdaError does not equal expected value");
0433         }
0434         if (muon.trk_phiError() != expectedMuonFloatingPointValues_[20] + offset) {
0435           throwWithMessage("analyzeMuons, trk_phiError does not equal expected value");
0436         }
0437         if (muon.trk_dsz() != expectedMuonFloatingPointValues_[21] + offset) {
0438           throwWithMessage("analyzeMuons, trk_dsz does not equal expected value");
0439         }
0440         if (muon.trk_dszError() != expectedMuonFloatingPointValues_[22] + offset) {
0441           throwWithMessage("analyzeMuons, trk_dszError does not equal expected value");
0442         }
0443         int j = 0;
0444         for (auto const& val : muon.vtxIndx()) {
0445           if (val != expectedMuonIntegralValues_[7] + iOffset + 10 * j) {
0446             throwWithMessage("analyzeMuons, vtxIndx does not contain expected value");
0447           }
0448           ++j;
0449         }
0450       }
0451       ++i;
0452     }
0453   }
0454 
0455   void TestReadRun2Scouting::analyzeParticles(edm::Event const& iEvent) const {
0456     auto const& particles = iEvent.get(particlesToken_);
0457     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0458     if (particles.size() != vectorSize) {
0459       throwWithMessage("analyzeParticles, particles does not have expected size");
0460     }
0461     unsigned int i = 0;
0462     for (auto const& particle : particles) {
0463       double offset = static_cast<double>(iEvent.id().event() + i);
0464       int iOffset = static_cast<int>(iEvent.id().event() + i);
0465 
0466       if (particle.pt() != expectedParticleFloatingPointValues_[0] + offset) {
0467         throwWithMessage("analyzeParticles, pt does not equal expected value");
0468       }
0469       if (particle.eta() != expectedParticleFloatingPointValues_[1] + offset) {
0470         throwWithMessage("analyzeParticles, eta does not equal expected value");
0471       }
0472       if (particle.phi() != expectedParticleFloatingPointValues_[2] + offset) {
0473         throwWithMessage("analyzeParticles, phi does not equal expected value");
0474       }
0475       if (particle.m() != expectedParticleFloatingPointValues_[3] + offset) {
0476         throwWithMessage("analyzeParticles, m does not equal expected value");
0477       }
0478       if (particle.pdgId() != expectedParticleIntegralValues_[0] + iOffset) {
0479         throwWithMessage("analyzeParticles, pdgId does not equal expected value");
0480       }
0481       if (particle.vertex() != expectedParticleIntegralValues_[1] + iOffset) {
0482         throwWithMessage("analyzeParticles, vertex does not equal expected value");
0483       }
0484       ++i;
0485     }
0486   }
0487 
0488   void TestReadRun2Scouting::analyzePFJets(edm::Event const& iEvent) const {
0489     auto const& pfJets = iEvent.get(pfJetsToken_);
0490     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0491     if (pfJets.size() != vectorSize) {
0492       throwWithMessage("analyzePFJets, pfJets does not have expected size");
0493     }
0494     unsigned int i = 0;
0495     for (auto const& pfJet : pfJets) {
0496       double offset = static_cast<double>(iEvent.id().event() + i);
0497       int iOffset = static_cast<int>(iEvent.id().event() + i);
0498 
0499       if (pfJet.pt() != expectedPFJetFloatingPointValues_[0] + offset) {
0500         throwWithMessage("analyzePFJets, pt does not equal expected value");
0501       }
0502       if (pfJet.eta() != expectedPFJetFloatingPointValues_[1] + offset) {
0503         throwWithMessage("analyzePFJets, eta does not equal expected value");
0504       }
0505       if (pfJet.phi() != expectedPFJetFloatingPointValues_[2] + offset) {
0506         throwWithMessage("analyzePFJets, phi does not equal expected value");
0507       }
0508       if (pfJet.m() != expectedPFJetFloatingPointValues_[3] + offset) {
0509         throwWithMessage("analyzePFJets, m does not equal expected value");
0510       }
0511       if (pfJet.jetArea() != expectedPFJetFloatingPointValues_[4] + offset) {
0512         throwWithMessage("analyzePFJets, jetArea does not equal expected value");
0513       }
0514       if (pfJet.chargedHadronEnergy() != expectedPFJetFloatingPointValues_[5] + offset) {
0515         throwWithMessage("analyzePFJets, chargedHadronEnergy does not equal expected value");
0516       }
0517       if (pfJet.neutralHadronEnergy() != expectedPFJetFloatingPointValues_[6] + offset) {
0518         throwWithMessage("analyzePFJets, neutralHadronEnergy does not equal expected value");
0519       }
0520       if (pfJet.photonEnergy() != expectedPFJetFloatingPointValues_[7] + offset) {
0521         throwWithMessage("analyzePFJets, photonEnergy does not equal expected value");
0522       }
0523       if (pfJet.electronEnergy() != expectedPFJetFloatingPointValues_[8] + offset) {
0524         throwWithMessage("analyzePFJets, electronEnergy does not equal expected value");
0525       }
0526       if (pfJet.muonEnergy() != expectedPFJetFloatingPointValues_[9] + offset) {
0527         throwWithMessage("analyzePFJets, muonEnergy does not equal expected value");
0528       }
0529       if (pfJet.HFHadronEnergy() != expectedPFJetFloatingPointValues_[10] + offset) {
0530         throwWithMessage("analyzePFJets, HFHadronEnergy does not equal expected value");
0531       }
0532       if (pfJet.HFEMEnergy() != expectedPFJetFloatingPointValues_[11] + offset) {
0533         throwWithMessage("analyzePFJets, HFEMEnergy does not equal expected value");
0534       }
0535       if (pfJet.chargedHadronMultiplicity() != expectedPFJetIntegralValues_[0] + iOffset) {
0536         throwWithMessage("analyzePFJets, chargedHadronMultiplicity does not equal expected value");
0537       }
0538       if (pfJet.neutralHadronMultiplicity() != expectedPFJetIntegralValues_[1] + iOffset) {
0539         throwWithMessage("analyzePFJets, neutralHadronMultiplicity does not equal expected value");
0540       }
0541       if (pfJet.photonMultiplicity() != expectedPFJetIntegralValues_[2] + iOffset) {
0542         throwWithMessage("analyzePFJets, photonMultiplicity does not equal expected value");
0543       }
0544       if (pfJet.electronMultiplicity() != expectedPFJetIntegralValues_[3] + iOffset) {
0545         throwWithMessage("analyzePFJets, electronMultiplicity does not equal expected value");
0546       }
0547       if (pfJet.muonMultiplicity() != expectedPFJetIntegralValues_[4] + iOffset) {
0548         throwWithMessage("analyzePFJets, muonMultiplicity does not equal expected value");
0549       }
0550       if (pfJet.HFHadronMultiplicity() != expectedPFJetIntegralValues_[5] + iOffset) {
0551         throwWithMessage("analyzePFJets, HFHadronMultiplicity does not equal expected value");
0552       }
0553       if (pfJet.HFEMMultiplicity() != expectedPFJetIntegralValues_[6] + iOffset) {
0554         throwWithMessage("analyzePFJets, HFEMMultiplicity does not equal expected value");
0555       }
0556       if (pfJet.HOEnergy() != expectedPFJetFloatingPointValues_[12] + offset) {
0557         throwWithMessage("analyzePFJets, HOEnergy does not equal expected value");
0558       }
0559       if (pfJet.csv() != expectedPFJetFloatingPointValues_[13] + offset) {
0560         throwWithMessage("analyzePFJets, csv does not equal expected value");
0561       }
0562       if (pfJet.mvaDiscriminator() != expectedPFJetFloatingPointValues_[14] + offset) {
0563         throwWithMessage("analyzePFJets, mvaDiscriminator does not equal expected value");
0564       }
0565       int j = 0;
0566       for (auto const& val : pfJet.constituents()) {
0567         if (val != expectedPFJetIntegralValues_[7] + iOffset + 10 * j) {
0568           throwWithMessage("analyzePFJets, constituents does not contain expected value");
0569         }
0570         ++j;
0571       }
0572       ++i;
0573     }
0574   }
0575 
0576   void TestReadRun2Scouting::analyzePhotons(edm::Event const& iEvent) const {
0577     auto const& photons = iEvent.get(photonsToken_);
0578     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0579     if (photons.size() != vectorSize) {
0580       throwWithMessage("analyzePhotons, photons does not have expected size");
0581     }
0582     unsigned int i = 0;
0583     for (auto const& photon : photons) {
0584       double offset = static_cast<double>(iEvent.id().event() + i);
0585 
0586       if (photon.pt() != expectedPhotonFloatingPointValues_[0] + offset) {
0587         throwWithMessage("analyzePhotons, pt does not equal expected value");
0588       }
0589       if (photon.eta() != expectedPhotonFloatingPointValues_[1] + offset) {
0590         throwWithMessage("analyzePhotons, eta does not equal expected value");
0591       }
0592       if (photon.phi() != expectedPhotonFloatingPointValues_[2] + offset) {
0593         throwWithMessage("analyzePhotons, phi does not equal expected value");
0594       }
0595       if (photon.m() != expectedPhotonFloatingPointValues_[3] + offset) {
0596         throwWithMessage("analyzePhotons, m does not equal expected value");
0597       }
0598       if (photon.sigmaIetaIeta() != expectedPhotonFloatingPointValues_[4] + offset) {
0599         throwWithMessage("analyzePhotons, sigmaIetaIeta does not equal expected value");
0600       }
0601       if (photon.hOverE() != expectedPhotonFloatingPointValues_[5] + offset) {
0602         throwWithMessage("analyzePhotons, hOverE does not equal expected value");
0603       }
0604       if (photon.ecalIso() != expectedPhotonFloatingPointValues_[6] + offset) {
0605         throwWithMessage("analyzePhotons, ecalIso does not equal expected value");
0606       }
0607       if (photon.hcalIso() != expectedPhotonFloatingPointValues_[7] + offset) {
0608         throwWithMessage("analyzePhotons, hcalIso does not equal expected value");
0609       }
0610       ++i;
0611     }
0612   }
0613 
0614   void TestReadRun2Scouting::analyzeTracks(edm::Event const& iEvent) const {
0615     if (inputTrackClassVersion_ < 2) {
0616       return;
0617     }
0618     auto const& tracks = iEvent.get(tracksToken_);
0619     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0620     if (tracks.size() != vectorSize) {
0621       throwWithMessage("analyzeTracks, tracks does not have expected size");
0622     }
0623     unsigned int i = 0;
0624     for (auto const& track : tracks) {
0625       double offset = static_cast<double>(iEvent.id().event() + i);
0626       int iOffset = static_cast<int>(iEvent.id().event() + i);
0627 
0628       if (track.tk_pt() != expectedTrackFloatingPointValues_[0] + offset) {
0629         throwWithMessage("analyzeTracks, tk_pt does not equal expected value");
0630       }
0631       if (track.tk_eta() != expectedTrackFloatingPointValues_[1] + offset) {
0632         throwWithMessage("analyzeTracks, tk_eta does not equal expected value");
0633       }
0634       if (track.tk_phi() != expectedTrackFloatingPointValues_[2] + offset) {
0635         throwWithMessage("analyzeTracks, tk_phi does not equal expected value");
0636       }
0637       if (track.tk_chi2() != expectedTrackFloatingPointValues_[3] + offset) {
0638         throwWithMessage("analyzeTracks, tk_chi2 does not equal expected value");
0639       }
0640       if (track.tk_ndof() != expectedTrackFloatingPointValues_[4] + offset) {
0641         throwWithMessage("analyzeTracks, tk_ndof does not equal expected value");
0642       }
0643       if (track.tk_charge() != expectedTrackIntegralValues_[0] + iOffset) {
0644         throwWithMessage("analyzeTracks, tk_charge does not equal expected value");
0645       }
0646       if (track.tk_dxy() != expectedTrackFloatingPointValues_[5] + offset) {
0647         throwWithMessage("analyzeTracks, tk_dxy does not equal expected value");
0648       }
0649       if (track.tk_dz() != expectedTrackFloatingPointValues_[6] + offset) {
0650         throwWithMessage("analyzeTracks, tk_dz does not equal expected value");
0651       }
0652       if (track.tk_nValidPixelHits() != expectedTrackIntegralValues_[1] + iOffset) {
0653         throwWithMessage("analyzeTracks, tk_nValidPixelHits does not equal expected value");
0654       }
0655       if (track.tk_nTrackerLayersWithMeasurement() != expectedTrackIntegralValues_[2] + iOffset) {
0656         throwWithMessage("analyzeTracks, tk_nTrackerLayersWithMeasurement does not equal expected value");
0657       }
0658       if (track.tk_nValidStripHits() != expectedTrackIntegralValues_[3] + iOffset) {
0659         throwWithMessage("analyzeTracks, tk_nValidStripHits does not equal expected value");
0660       }
0661       if (track.tk_qoverp() != expectedTrackFloatingPointValues_[7] + offset) {
0662         throwWithMessage("analyzeTracks, tk_qoverp does not equal expected value");
0663       }
0664       if (track.tk_lambda() != expectedTrackFloatingPointValues_[8] + offset) {
0665         throwWithMessage("analyzeTracks, tk_lambda does not equal expected value");
0666       }
0667       if (track.tk_dxy_Error() != expectedTrackFloatingPointValues_[9] + offset) {
0668         throwWithMessage("analyzeTracks, tk_dxy_Error does not equal expected value");
0669       }
0670       if (track.tk_dz_Error() != expectedTrackFloatingPointValues_[10] + offset) {
0671         throwWithMessage("analyzeTracks, tk_dz_Error does not equal expected value");
0672       }
0673       if (track.tk_qoverp_Error() != expectedTrackFloatingPointValues_[11] + offset) {
0674         throwWithMessage("analyzeTracks, tk_qoverp_Error does not equal expected value");
0675       }
0676       if (track.tk_lambda_Error() != expectedTrackFloatingPointValues_[12] + offset) {
0677         throwWithMessage("analyzeTracks, tk_lambda_Error does not equal expected value");
0678       }
0679       if (track.tk_phi_Error() != expectedTrackFloatingPointValues_[13] + offset) {
0680         throwWithMessage("analyzeTracks, tk_phi_Error does not equal expected value");
0681       }
0682       if (track.tk_dsz() != expectedTrackFloatingPointValues_[14] + offset) {
0683         throwWithMessage("analyzeTracks, tk_dsz does not equal expected value");
0684       }
0685       if (track.tk_dsz_Error() != expectedTrackFloatingPointValues_[15] + offset) {
0686         throwWithMessage("analyzeTracks, tk_dsz_Error does not equal expected value");
0687       }
0688       ++i;
0689     }
0690   }
0691 
0692   void TestReadRun2Scouting::analyzeVertexes(edm::Event const& iEvent) const {
0693     auto const& vertexes = iEvent.get(vertexesToken_);
0694     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0695     if (vertexes.size() != vectorSize) {
0696       throwWithMessage("analyzeVertexes, vertexes does not have expected size");
0697     }
0698     unsigned int i = 0;
0699     for (auto const& vertex : vertexes) {
0700       double offset = static_cast<double>(iEvent.id().event() + i);
0701       int iOffset = static_cast<int>(iEvent.id().event() + i);
0702 
0703       if (vertex.x() != expectedVertexFloatingPointValues_[0] + offset) {
0704         throwWithMessage("analyzeVertexes, x does not equal expected value");
0705       }
0706       if (vertex.y() != expectedVertexFloatingPointValues_[1] + offset) {
0707         throwWithMessage("analyzeVertexes, y does not equal expected value");
0708       }
0709       if (vertex.z() != expectedVertexFloatingPointValues_[2] + offset) {
0710         throwWithMessage("analyzeVertexes, z does not equal expected value");
0711       }
0712       if (vertex.zError() != expectedVertexFloatingPointValues_[3] + offset) {
0713         throwWithMessage("analyzeVertexes, zError does not equal expected value");
0714       }
0715       if (inputVertexClassVersion_ > 2) {
0716         if (vertex.xError() != expectedVertexFloatingPointValues_[4] + offset) {
0717           throwWithMessage("analyzeVertexes, xError does not equal expected value");
0718         }
0719         if (vertex.yError() != expectedVertexFloatingPointValues_[5] + offset) {
0720           throwWithMessage("analyzeVertexes, yError does not equal expected value");
0721         }
0722         if (vertex.tracksSize() != expectedVertexIntegralValues_[0] + iOffset) {
0723           throwWithMessage("analyzeVertexes, tracksSize does not equal expected value");
0724         }
0725         if (vertex.chi2() != expectedVertexFloatingPointValues_[6] + offset) {
0726           throwWithMessage("analyzeVertexes, chi2 does not equal expected value");
0727         }
0728         if (vertex.ndof() != expectedVertexIntegralValues_[1] + iOffset) {
0729           throwWithMessage("analyzeVertexes, ndof does not equal expected value");
0730         }
0731         if (vertex.isValidVtx() != static_cast<bool>((expectedVertexIntegralValues_[2] + iOffset) % 2)) {
0732           throwWithMessage("analyzeVertexes, isValidVtx does not equal expected value");
0733         }
0734       }
0735       ++i;
0736     }
0737   }
0738 
0739   void TestReadRun2Scouting::throwWithMessageFromConstructor(const char* msg) const {
0740     throw cms::Exception("TestFailure") << "TestReadRun2Scouting constructor, " << msg;
0741   }
0742 
0743   void TestReadRun2Scouting::throwWithMessage(const char* msg) const {
0744     throw cms::Exception("TestFailure") << "TestReadRun2Scouting::analyze, " << msg;
0745   }
0746 
0747 }  // namespace edmtest
0748 
0749 using edmtest::TestReadRun2Scouting;
0750 DEFINE_FWK_MODULE(TestReadRun2Scouting);