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:      TestReadRun3Scouting
0005 //
0006 /**\class edmtest::TestReadRun3Scouting
0007   Description: Used as part of tests that ensure the run 3 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 3 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:  18 May 2023
0017 
0018 #include "DataFormats/Scouting/interface/Run3ScoutingCaloJet.h"
0019 #include "DataFormats/Scouting/interface/Run3ScoutingElectron.h"
0020 #include "DataFormats/Scouting/interface/Run3ScoutingHitPatternPOD.h"
0021 #include "DataFormats/Scouting/interface/Run3ScoutingMuon.h"
0022 #include "DataFormats/Scouting/interface/Run3ScoutingParticle.h"
0023 #include "DataFormats/Scouting/interface/Run3ScoutingPFJet.h"
0024 #include "DataFormats/Scouting/interface/Run3ScoutingPhoton.h"
0025 #include "DataFormats/Scouting/interface/Run3ScoutingTrack.h"
0026 #include "DataFormats/Scouting/interface/Run3ScoutingVertex.h"
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/Frameworkfwd.h"
0029 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0034 #include "FWCore/Utilities/interface/EDGetToken.h"
0035 #include "FWCore/Utilities/interface/Exception.h"
0036 #include "FWCore/Utilities/interface/InputTag.h"
0037 #include "FWCore/Utilities/interface/StreamID.h"
0038 
0039 #include <vector>
0040 
0041 namespace edmtest {
0042 
0043   class TestReadRun3Scouting : public edm::global::EDAnalyzer<> {
0044   public:
0045     TestReadRun3Scouting(edm::ParameterSet const&);
0046     void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;
0047     static void fillDescriptions(edm::ConfigurationDescriptions&);
0048 
0049   private:
0050     void analyzeCaloJets(edm::Event const&) const;
0051     void analyzeElectrons(edm::Event const&) const;
0052     void analyzeMuons(edm::Event const&) const;
0053     void analyzeParticles(edm::Event const&) const;
0054     void analyzePFJets(edm::Event const&) const;
0055     void analyzePhotons(edm::Event const&) const;
0056     void analyzeTracks(edm::Event const&) const;
0057     void analyzeVertexes(edm::Event const&) const;
0058 
0059     void throwWithMessageFromConstructor(const char*) const;
0060     void throwWithMessage(const char*) const;
0061 
0062     // These expected values are meaningless other than we use them
0063     // to check that values read from persistent storage match the values
0064     // we know were written.
0065 
0066     const std::vector<double> expectedCaloJetsValues_;
0067     const edm::EDGetTokenT<std::vector<Run3ScoutingCaloJet>> caloJetsToken_;
0068 
0069     const int inputElectronClassVersion_;
0070     const std::vector<double> expectedElectronFloatingPointValues_;
0071     const std::vector<int> expectedElectronIntegralValues_;
0072     const edm::EDGetTokenT<std::vector<Run3ScoutingElectron>> electronsToken_;
0073 
0074     const std::vector<double> expectedMuonFloatingPointValues_;
0075     const std::vector<int> expectedMuonIntegralValues_;
0076     const edm::EDGetTokenT<std::vector<Run3ScoutingMuon>> muonsToken_;
0077 
0078     const std::vector<double> expectedParticleFloatingPointValues_;
0079     const std::vector<int> expectedParticleIntegralValues_;
0080     const edm::EDGetTokenT<std::vector<Run3ScoutingParticle>> particlesToken_;
0081 
0082     const std::vector<double> expectedPFJetFloatingPointValues_;
0083     const std::vector<int> expectedPFJetIntegralValues_;
0084     const edm::EDGetTokenT<std::vector<Run3ScoutingPFJet>> pfJetsToken_;
0085 
0086     const int inputPhotonClassVersion_;
0087     const std::vector<double> expectedPhotonFloatingPointValues_;
0088     const std::vector<int> expectedPhotonIntegralValues_;
0089     const edm::EDGetTokenT<std::vector<Run3ScoutingPhoton>> photonsToken_;
0090 
0091     const std::vector<double> expectedTrackFloatingPointValues_;
0092     const std::vector<int> expectedTrackIntegralValues_;
0093     const edm::EDGetTokenT<std::vector<Run3ScoutingTrack>> tracksToken_;
0094 
0095     const int inputVertexClassVersion_;
0096     const std::vector<double> expectedVertexFloatingPointValues_;
0097     const std::vector<int> expectedVertexIntegralValues_;
0098     const edm::EDGetTokenT<std::vector<Run3ScoutingVertex>> vertexesToken_;
0099   };
0100 
0101   TestReadRun3Scouting::TestReadRun3Scouting(edm::ParameterSet const& iPSet)
0102       : expectedCaloJetsValues_(iPSet.getParameter<std::vector<double>>("expectedCaloJetsValues")),
0103         caloJetsToken_(consumes(iPSet.getParameter<edm::InputTag>("caloJetsTag"))),
0104         inputElectronClassVersion_(iPSet.getParameter<int>("electronClassVersion")),
0105         expectedElectronFloatingPointValues_(
0106             iPSet.getParameter<std::vector<double>>("expectedElectronFloatingPointValues")),
0107         expectedElectronIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedElectronIntegralValues")),
0108         electronsToken_(consumes(iPSet.getParameter<edm::InputTag>("electronsTag"))),
0109         expectedMuonFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedMuonFloatingPointValues")),
0110         expectedMuonIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedMuonIntegralValues")),
0111         muonsToken_(consumes(iPSet.getParameter<edm::InputTag>("muonsTag"))),
0112         expectedParticleFloatingPointValues_(
0113             iPSet.getParameter<std::vector<double>>("expectedParticleFloatingPointValues")),
0114         expectedParticleIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedParticleIntegralValues")),
0115         particlesToken_(consumes(iPSet.getParameter<edm::InputTag>("particlesTag"))),
0116         expectedPFJetFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedPFJetFloatingPointValues")),
0117         expectedPFJetIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedPFJetIntegralValues")),
0118         pfJetsToken_(consumes(iPSet.getParameter<edm::InputTag>("pfJetsTag"))),
0119         inputPhotonClassVersion_(iPSet.getParameter<int>("photonClassVersion")),
0120         expectedPhotonFloatingPointValues_(
0121             iPSet.getParameter<std::vector<double>>("expectedPhotonFloatingPointValues")),
0122         expectedPhotonIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedPhotonIntegralValues")),
0123         photonsToken_(consumes(iPSet.getParameter<edm::InputTag>("photonsTag"))),
0124         expectedTrackFloatingPointValues_(iPSet.getParameter<std::vector<double>>("expectedTrackFloatingPointValues")),
0125         expectedTrackIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedTrackIntegralValues")),
0126         tracksToken_(consumes(iPSet.getParameter<edm::InputTag>("tracksTag"))),
0127         inputVertexClassVersion_(iPSet.getParameter<int>("vertexClassVersion")),
0128         expectedVertexFloatingPointValues_(
0129             iPSet.getParameter<std::vector<double>>("expectedVertexFloatingPointValues")),
0130         expectedVertexIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedVertexIntegralValues")),
0131         vertexesToken_(consumes(iPSet.getParameter<edm::InputTag>("vertexesTag"))) {
0132     if (expectedCaloJetsValues_.size() != 16) {
0133       throwWithMessageFromConstructor("test configuration error, expectedCaloJetsValues must have size 16");
0134     }
0135     if (expectedElectronFloatingPointValues_.size() != 33) {
0136       throwWithMessageFromConstructor(
0137           "test configuration error, expectedElectronFloatingPointValues must have size 33");
0138     }
0139     if (expectedElectronIntegralValues_.size() != 8) {
0140       throwWithMessageFromConstructor("test configuration error, expectedElectronIntegralValues must have size 8");
0141     }
0142     if (expectedMuonFloatingPointValues_.size() != 37) {
0143       throwWithMessageFromConstructor("test configuration error, expectedMuonFloatingPointValues must have size 37");
0144     }
0145     if (expectedMuonIntegralValues_.size() != 26) {
0146       throwWithMessageFromConstructor("test configuration error, expectedMuonIntegralValues must have size 26");
0147     }
0148     if (expectedParticleFloatingPointValues_.size() != 11) {
0149       throwWithMessageFromConstructor(
0150           "test configuration error, expectedParticleFloatingPointValues must have size 11");
0151     }
0152     if (expectedParticleIntegralValues_.size() != 5) {
0153       throwWithMessageFromConstructor("test configuration error, expectedParticleIntegralValues must have size 5");
0154     }
0155     if (expectedPFJetFloatingPointValues_.size() != 15) {
0156       throwWithMessageFromConstructor("test configuration error, expectedPFJetFloatingPointValues must have size 15");
0157     }
0158     if (expectedPFJetIntegralValues_.size() != 8) {
0159       throwWithMessageFromConstructor("test configuration error, expectedPFJetIntegralValues must have size 8");
0160     }
0161     if (expectedPhotonFloatingPointValues_.size() != 17) {
0162       throwWithMessageFromConstructor("test configuration error, expectedPhotonFloatingPointValues must have size 17");
0163     }
0164     if (expectedPhotonIntegralValues_.size() != 5) {
0165       throwWithMessageFromConstructor("test configuration error, expectedPhotonIntegralValues must have size 5");
0166     }
0167     if (expectedTrackFloatingPointValues_.size() != 29) {
0168       throwWithMessageFromConstructor("test configuration error, expectedTrackFloatingPointValues must have size 29");
0169     }
0170     if (expectedTrackIntegralValues_.size() != 5) {
0171       throwWithMessageFromConstructor("test configuration error, expectedTrackIntegralValues must have size 5");
0172     }
0173     if (expectedVertexFloatingPointValues_.size() != 10) {
0174       throwWithMessageFromConstructor("test configuration error, expectedVertexFloatingPointValues must have size 10");
0175     }
0176     if (expectedVertexIntegralValues_.size() != 3) {
0177       throwWithMessageFromConstructor("test configuration error, expectedVertexIntegralValues must have size 3");
0178     }
0179   }
0180 
0181   void TestReadRun3Scouting::analyze(edm::StreamID, edm::Event const& iEvent, edm::EventSetup const&) const {
0182     analyzeCaloJets(iEvent);
0183     analyzeElectrons(iEvent);
0184     analyzeMuons(iEvent);
0185     analyzeParticles(iEvent);
0186     analyzePFJets(iEvent);
0187     analyzePhotons(iEvent);
0188     analyzeTracks(iEvent);
0189     analyzeVertexes(iEvent);
0190   }
0191 
0192   void TestReadRun3Scouting::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0193     edm::ParameterSetDescription desc;
0194     desc.add<std::vector<double>>("expectedCaloJetsValues");
0195     desc.add<edm::InputTag>("caloJetsTag");
0196     desc.add<int>("electronClassVersion");
0197     desc.add<std::vector<double>>("expectedElectronFloatingPointValues");
0198     desc.add<std::vector<int>>("expectedElectronIntegralValues");
0199     desc.add<edm::InputTag>("electronsTag");
0200     desc.add<std::vector<double>>("expectedMuonFloatingPointValues");
0201     desc.add<std::vector<int>>("expectedMuonIntegralValues");
0202     desc.add<edm::InputTag>("muonsTag");
0203     desc.add<std::vector<double>>("expectedParticleFloatingPointValues");
0204     desc.add<std::vector<int>>("expectedParticleIntegralValues");
0205     desc.add<edm::InputTag>("particlesTag");
0206     desc.add<std::vector<double>>("expectedPFJetFloatingPointValues");
0207     desc.add<std::vector<int>>("expectedPFJetIntegralValues");
0208     desc.add<edm::InputTag>("pfJetsTag");
0209     desc.add<int>("photonClassVersion");
0210     desc.add<std::vector<double>>("expectedPhotonFloatingPointValues");
0211     desc.add<std::vector<int>>("expectedPhotonIntegralValues");
0212     desc.add<edm::InputTag>("photonsTag");
0213     desc.add<std::vector<double>>("expectedTrackFloatingPointValues");
0214     desc.add<std::vector<int>>("expectedTrackIntegralValues");
0215     desc.add<edm::InputTag>("tracksTag");
0216     desc.add<int>("vertexClassVersion");
0217     desc.add<std::vector<double>>("expectedVertexFloatingPointValues");
0218     desc.add<std::vector<int>>("expectedVertexIntegralValues");
0219     desc.add<edm::InputTag>("vertexesTag");
0220     descriptions.addDefault(desc);
0221   }
0222 
0223   void TestReadRun3Scouting::analyzeCaloJets(edm::Event const& iEvent) const {
0224     auto const& caloJets = iEvent.get(caloJetsToken_);
0225     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0226     if (caloJets.size() != vectorSize) {
0227       throwWithMessage("analyzeCaloJets, caloJets does not have expected size");
0228     }
0229     unsigned int i = 0;
0230     for (auto const& caloJet : caloJets) {
0231       double offset = static_cast<double>(iEvent.id().event() + i);
0232 
0233       if (caloJet.pt() != expectedCaloJetsValues_[0] + offset) {
0234         throwWithMessage("analyzeCaloJets, pt does not equal expected value");
0235       }
0236       if (caloJet.eta() != expectedCaloJetsValues_[1] + offset) {
0237         throwWithMessage("analyzeCaloJets, eta does not equal expected value");
0238       }
0239       if (caloJet.phi() != expectedCaloJetsValues_[2] + offset) {
0240         throwWithMessage("analyzeCaloJets, phi does not equal expected value");
0241       }
0242       if (caloJet.m() != expectedCaloJetsValues_[3] + offset) {
0243         throwWithMessage("analyzeCaloJets, m does not equal expected value");
0244       }
0245       if (caloJet.jetArea() != expectedCaloJetsValues_[4] + offset) {
0246         throwWithMessage("analyzeCaloJets, jetArea does not equal expected value");
0247       }
0248       if (caloJet.maxEInEmTowers() != expectedCaloJetsValues_[5] + offset) {
0249         throwWithMessage("analyzeCaloJets,  maxEInEmTowers() does not equal expected value");
0250       }
0251       if (caloJet.maxEInHadTowers() != expectedCaloJetsValues_[6] + offset) {
0252         throwWithMessage("analyzeCaloJets,  maxEInHadTowers does not equal expected value");
0253       }
0254       if (caloJet.hadEnergyInHB() != expectedCaloJetsValues_[7] + offset) {
0255         throwWithMessage("analyzeCaloJets, hadEnergyInHB does not equal expected value");
0256       }
0257       if (caloJet.hadEnergyInHE() != expectedCaloJetsValues_[8] + offset) {
0258         throwWithMessage("analyzeCaloJets, hadEnergyInHE does not equal expected value");
0259       }
0260       if (caloJet.hadEnergyInHF() != expectedCaloJetsValues_[9] + offset) {
0261         throwWithMessage("analyzeCaloJets, hadEnergyInHF does not equal expected value");
0262       }
0263       if (caloJet.emEnergyInEB() != expectedCaloJetsValues_[10] + offset) {
0264         throwWithMessage("analyzeCaloJets, emEnergyInEB does not equal expected value");
0265       }
0266       if (caloJet.emEnergyInEE() != expectedCaloJetsValues_[11] + offset) {
0267         throwWithMessage("analyzeCaloJets, emEnergyInEE does not equal expected value");
0268       }
0269       if (caloJet.emEnergyInHF() != expectedCaloJetsValues_[12] + offset) {
0270         throwWithMessage("analyzeCaloJets, emEnergyInHF does not equal expected value");
0271       }
0272       if (caloJet.towersArea() != expectedCaloJetsValues_[13] + offset) {
0273         throwWithMessage("analyzeCaloJets, towersArea does not equal expected value");
0274       }
0275       if (caloJet.mvaDiscriminator() != expectedCaloJetsValues_[14] + offset) {
0276         throwWithMessage("analyzeCaloJets,  mvaDiscriminator does not equal expected value");
0277       }
0278       if (caloJet.btagDiscriminator() != expectedCaloJetsValues_[15] + offset) {
0279         throwWithMessage("analyzeCaloJets,  btagDiscriminator does not equal expected value");
0280       }
0281       ++i;
0282     }
0283   }
0284 
0285   void TestReadRun3Scouting::analyzeElectrons(edm::Event const& iEvent) const {
0286     auto const& electrons = iEvent.get(electronsToken_);
0287     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0288     if (electrons.size() != vectorSize) {
0289       throwWithMessage("analyzeElectrons, electrons does not have expected size");
0290     }
0291     unsigned int i = 0;
0292     for (auto const& electron : electrons) {
0293       double offset = static_cast<double>(iEvent.id().event() + i);
0294       int iOffset = static_cast<int>(iEvent.id().event() + i);
0295 
0296       if (electron.pt() != expectedElectronFloatingPointValues_[0] + offset) {
0297         throwWithMessage("analyzeElectrons, pt does not equal expected value");
0298       }
0299       if (electron.eta() != expectedElectronFloatingPointValues_[1] + offset) {
0300         throwWithMessage("analyzeElectrons, eta does not equal expected value");
0301       }
0302       if (electron.phi() != expectedElectronFloatingPointValues_[2] + offset) {
0303         throwWithMessage("analyzeElectrons, phi does not equal expected value");
0304       }
0305       if (electron.m() != expectedElectronFloatingPointValues_[3] + offset) {
0306         throwWithMessage("analyzeElectrons, m does not equal expected value");
0307       }
0308       if (inputElectronClassVersion_ == 5) {
0309         if (electron.trkd0().size() != 1) {
0310           throwWithMessage("analyzeElectrons, trkd0 does not have expected size");
0311         }
0312         if (electron.trkd0()[0] != expectedElectronFloatingPointValues_[4] + offset) {
0313           throwWithMessage("analyzeElectrons, d0 does not equal expected value");
0314         }
0315         if (electron.trkdz().size() != 1) {
0316           throwWithMessage("analyzeElectrons, trkdz does not have expected size");
0317         }
0318         if (electron.trkdz()[0] != expectedElectronFloatingPointValues_[5] + offset) {
0319           throwWithMessage("analyzeElectrons,  dz does not equal expected value");
0320         }
0321       }
0322       if (electron.dEtaIn() != expectedElectronFloatingPointValues_[6] + offset) {
0323         throwWithMessage("analyzeElectrons,  dEtaIn does not equal expected value");
0324       }
0325       if (electron.dPhiIn() != expectedElectronFloatingPointValues_[7] + offset) {
0326         throwWithMessage("analyzeElectrons, dPhiIn does not equal expected value");
0327       }
0328       if (electron.sigmaIetaIeta() != expectedElectronFloatingPointValues_[8] + offset) {
0329         throwWithMessage("analyzeElectrons, sigmaIetaIeta does not equal expected value");
0330       }
0331       if (electron.hOverE() != expectedElectronFloatingPointValues_[9] + offset) {
0332         throwWithMessage("analyzeElectrons, hOverE does not equal expected value");
0333       }
0334       if (electron.ooEMOop() != expectedElectronFloatingPointValues_[10] + offset) {
0335         throwWithMessage("analyzeElectrons, ooEMOop does not equal expected value");
0336       }
0337       if (electron.missingHits() != expectedElectronIntegralValues_[0] + iOffset) {
0338         throwWithMessage("analyzeElectrons, missingHits does not equal expected value");
0339       }
0340       if (inputElectronClassVersion_ == 5) {
0341         if (electron.trkcharge().size() != 1) {
0342           throwWithMessage("analyzeElectrons, trkcharge does not have expected size");
0343         }
0344         if (electron.trkcharge()[0] != expectedElectronIntegralValues_[1] + iOffset) {
0345           throwWithMessage("analyzeElectrons, charge does not equal expected value");
0346         }
0347       }
0348       if (electron.ecalIso() != expectedElectronFloatingPointValues_[11] + offset) {
0349         throwWithMessage("analyzeElectrons, ecalIso does not equal expected value");
0350       }
0351       if (electron.hcalIso() != expectedElectronFloatingPointValues_[12] + offset) {
0352         throwWithMessage("analyzeElectrons, hcalIso does not equal expected value");
0353       }
0354       if (electron.trackIso() != expectedElectronFloatingPointValues_[13] + offset) {
0355         throwWithMessage("analyzeElectrons, trackIso does not equal expected value");
0356       }
0357       if (electron.r9() != expectedElectronFloatingPointValues_[14] + offset) {
0358         throwWithMessage("analyzeElectrons, r9 does not equal expected value");
0359       }
0360       if (electron.sMin() != expectedElectronFloatingPointValues_[15] + offset) {
0361         throwWithMessage("analyzeElectrons, sMin does not equal expected value");
0362       }
0363       if (electron.sMaj() != expectedElectronFloatingPointValues_[16] + offset) {
0364         throwWithMessage("analyzeElectrons, sMaj does not equal expected value");
0365       }
0366       if (electron.seedId() != static_cast<unsigned int>(expectedElectronIntegralValues_[2] + iOffset)) {
0367         throwWithMessage("analyzeElectrons, seedId does not equal expected value");
0368       }
0369       if (electron.energyMatrix().size() != vectorSize) {
0370         throwWithMessage("analyzeElectrons, energyMatrix does not have expected size");
0371       }
0372       unsigned int j = 0;
0373       for (auto const& val : electron.energyMatrix()) {
0374         if (val != expectedElectronFloatingPointValues_[17] + offset + 10 * j) {
0375           throwWithMessage("analyzeElectrons, energyMatrix does not contain expected value");
0376         }
0377         ++j;
0378       }
0379       if (electron.detIds().size() != vectorSize) {
0380         throwWithMessage("analyzeElectrons, detIds does not have expected size");
0381       }
0382       j = 0;
0383       for (auto const& val : electron.detIds()) {
0384         if (val != expectedElectronIntegralValues_[3] + iOffset + 10 * j) {
0385           throwWithMessage("analyzeElectrons, detIds does not contain expected value");
0386         }
0387         ++j;
0388       }
0389       if (electron.timingMatrix().size() != vectorSize) {
0390         throwWithMessage("analyzeElectrons, timingMatrix does not have expected size");
0391       }
0392       j = 0;
0393       for (auto const& val : electron.timingMatrix()) {
0394         if (val != expectedElectronFloatingPointValues_[18] + offset + 10 * j) {
0395           throwWithMessage("analyzeElectrons, timingMatrix does not contain expected value");
0396         }
0397         ++j;
0398       }
0399       if (electron.rechitZeroSuppression() != static_cast<bool>((expectedElectronIntegralValues_[4] + iOffset) % 2)) {
0400         throwWithMessage("analyzeElectrons, rechitZeroSuppression does not equal expected value");
0401       }
0402       if (inputElectronClassVersion_ == 6 || inputElectronClassVersion_ == 7) {
0403         if (electron.trkd0().size() != vectorSize) {
0404           throwWithMessage("analyzeElectrons, trkd0 does not have expected size");
0405         }
0406         j = 0;
0407         for (auto const& val : electron.trkd0()) {
0408           if (val != expectedElectronFloatingPointValues_[19] + offset + 10 * j) {
0409             throwWithMessage("analyzeElectrons, trkd0 does not contain expected value");
0410           }
0411           ++j;
0412         }
0413         if (electron.trkdz().size() != vectorSize) {
0414           throwWithMessage("analyzeElectrons, trkdz does not have expected size");
0415         }
0416         j = 0;
0417         for (auto const& val : electron.trkdz()) {
0418           if (val != expectedElectronFloatingPointValues_[20] + offset + 10 * j) {
0419             throwWithMessage("analyzeElectrons, trkdz does not contain expected value");
0420           }
0421           ++j;
0422         }
0423         if (electron.trkpt().size() != vectorSize) {
0424           throwWithMessage("analyzeElectrons, trkpt does not have expected size");
0425         }
0426         j = 0;
0427         for (auto const& val : electron.trkpt()) {
0428           if (val != expectedElectronFloatingPointValues_[21] + offset + 10 * j) {
0429             throwWithMessage("analyzeElectrons, trkpt does not contain expected value");
0430           }
0431           ++j;
0432         }
0433         if (electron.trketa().size() != vectorSize) {
0434           throwWithMessage("analyzeElectrons, trketa does not have expected size");
0435         }
0436         j = 0;
0437         for (auto const& val : electron.trketa()) {
0438           if (val != expectedElectronFloatingPointValues_[22] + offset + 10 * j) {
0439             throwWithMessage("analyzeElectrons, trketa does not contain expected value");
0440           }
0441           ++j;
0442         }
0443         if (electron.trkphi().size() != vectorSize) {
0444           throwWithMessage("analyzeElectrons, trkphi does not have expected size");
0445         }
0446         j = 0;
0447         for (auto const& val : electron.trkphi()) {
0448           if (val != expectedElectronFloatingPointValues_[23] + offset + 10 * j) {
0449             throwWithMessage("analyzeElectrons, trkphi does not contain expected value");
0450           }
0451           ++j;
0452         }
0453         if (electron.trkchi2overndf().size() != vectorSize) {
0454           throwWithMessage("analyzeElectrons, trkchi2overndf does not have expected size");
0455         }
0456         j = 0;
0457         for (auto const& val : electron.trkchi2overndf()) {
0458           if (val != expectedElectronFloatingPointValues_[24] + offset + 10 * j) {
0459             throwWithMessage("analyzeElectrons, trkchi2overndf does not contain expected value");
0460           }
0461           ++j;
0462         }
0463         if (electron.trkcharge().size() != vectorSize) {
0464           throwWithMessage("analyzeElectrons, trkcharge does not have expected size");
0465         }
0466         j = 0;
0467         for (auto const& val : electron.trkcharge()) {
0468           if (val != static_cast<int>(expectedElectronIntegralValues_[5] + iOffset + 10 * j)) {
0469             throwWithMessage("analyzeElectrons, trkcharge does not contain expected value");
0470           }
0471           ++j;
0472         }
0473       }
0474       if (inputElectronClassVersion_ == 7) {
0475         if (electron.rawEnergy() != expectedElectronFloatingPointValues_[25] + offset) {
0476           throwWithMessage("analyzeElectrons, rawEnergy does not equal expected value");
0477         }
0478         if (electron.preshowerEnergy() != expectedElectronFloatingPointValues_[26] + offset) {
0479           throwWithMessage("analyzeElectrons, preshowerEnergy does not equal expected value");
0480         }
0481         if (electron.corrEcalEnergyError() != expectedElectronFloatingPointValues_[27] + offset) {
0482           throwWithMessage("analyzeElectrons, corrEcalEnergyError does not equal expected value");
0483         }
0484         if (electron.trkpMode().size() != vectorSize) {
0485           throwWithMessage("analyzeElectrons, trkpMode does not have expected size");
0486         }
0487         j = 0;
0488         for (auto const& val : electron.trkpMode()) {
0489           if (val != expectedElectronFloatingPointValues_[28] + offset + 10 * j) {
0490             throwWithMessage("analyzeElectrons, trkpMode does not contain expected value");
0491           }
0492           ++j;
0493         }
0494         if (electron.trketaMode().size() != vectorSize) {
0495           throwWithMessage("analyzeElectrons, trketaMode does not have expected size");
0496         }
0497         j = 0;
0498         for (auto const& val : electron.trketaMode()) {
0499           if (val != expectedElectronFloatingPointValues_[29] + offset + 10 * j) {
0500             throwWithMessage("analyzeElectrons, trketaMode does not contain expected value");
0501           }
0502           ++j;
0503         }
0504         if (electron.trkphiMode().size() != vectorSize) {
0505           throwWithMessage("analyzeElectrons, trkphiMode does not have expected size");
0506         }
0507         j = 0;
0508         for (auto const& val : electron.trkphiMode()) {
0509           if (val != expectedElectronFloatingPointValues_[30] + offset + 10 * j) {
0510             throwWithMessage("analyzeElectrons, trkphiMode does not contain expected value");
0511           }
0512           ++j;
0513         }
0514         if (electron.trkqoverpModeError().size() != vectorSize) {
0515           throwWithMessage("analyzeElectrons, trkqoverpModeError does not have expected size");
0516         }
0517         j = 0;
0518         for (auto const& val : electron.trkqoverpModeError()) {
0519           if (val != expectedElectronFloatingPointValues_[31] + offset + 10 * j) {
0520             throwWithMessage("analyzeElectrons, trkqoverpModeError does not contain expected value");
0521           }
0522           ++j;
0523         }
0524         if (electron.trackfbrem() != expectedElectronFloatingPointValues_[32] + offset) {
0525           throwWithMessage("analyzeElectrons, trackfbrem does not equal expected value");
0526         }
0527         if (electron.nClusters() != static_cast<unsigned int>(expectedElectronIntegralValues_[6] + iOffset)) {
0528           throwWithMessage("analyzeElectrons, nClusters does not equal expected value");
0529         }
0530         if (electron.nCrystals() != static_cast<unsigned int>(expectedElectronIntegralValues_[7] + iOffset)) {
0531           throwWithMessage("analyzeElectrons, nCrystals does not equal expected value");
0532         }
0533       }
0534       ++i;
0535     }
0536   }
0537 
0538   void TestReadRun3Scouting::analyzeMuons(edm::Event const& iEvent) const {
0539     auto const& muons = iEvent.get(muonsToken_);
0540     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0541     if (muons.size() != vectorSize) {
0542       throwWithMessage("analyzeMuons, muons does not have expected size");
0543     }
0544     unsigned int i = 0;
0545     for (auto const& muon : muons) {
0546       double offset = static_cast<double>(iEvent.id().event() + i);
0547       int iOffset = static_cast<int>(iEvent.id().event() + i);
0548 
0549       if (muon.pt() != expectedMuonFloatingPointValues_[0] + offset) {
0550         throwWithMessage("analyzeMuons, pt does not equal expected value");
0551       }
0552       if (muon.eta() != expectedMuonFloatingPointValues_[1] + offset) {
0553         throwWithMessage("analyzeMuons, eta does not equal expected value");
0554       }
0555       if (muon.phi() != expectedMuonFloatingPointValues_[2] + offset) {
0556         throwWithMessage("analyzeMuons, phi does not equal expected value");
0557       }
0558       if (muon.m() != expectedMuonFloatingPointValues_[3] + offset) {
0559         throwWithMessage("analyzeMuons, m does not equal expected value");
0560       }
0561       if (muon.type() != static_cast<unsigned int>(expectedMuonIntegralValues_[0] + iOffset)) {
0562         throwWithMessage("analyzeMuons, type does not equal expected value");
0563       }
0564       if (muon.charge() != expectedMuonIntegralValues_[1] + iOffset) {
0565         throwWithMessage("analyzeMuons, charge does not equal expected value");
0566       }
0567       if (muon.normalizedChi2() != expectedMuonFloatingPointValues_[4] + offset) {
0568         throwWithMessage("analyzeMuons,  normalizedChi2 does not equal expected value");
0569       }
0570       if (muon.ecalIso() != expectedMuonFloatingPointValues_[5] + offset) {
0571         throwWithMessage("analyzeMuons, ecalIso does not equal expected value");
0572       }
0573       if (muon.hcalIso() != expectedMuonFloatingPointValues_[6] + offset) {
0574         throwWithMessage("analyzeMuons, hcalIso does not equal expected value");
0575       }
0576       if (muon.trackIso() != expectedMuonFloatingPointValues_[7] + offset) {
0577         throwWithMessage("analyzeMuons, trackIso does not equal expected value");
0578       }
0579       if (muon.nValidStandAloneMuonHits() != expectedMuonIntegralValues_[2] + iOffset) {
0580         throwWithMessage("analyzeMuons, nValidStandAloneMuonHits does not equal expected value");
0581       }
0582       if (muon.nStandAloneMuonMatchedStations() != expectedMuonIntegralValues_[3] + iOffset) {
0583         throwWithMessage("analyzeMuons, nStandAloneMuonMatchedStations does not equal expected value");
0584       }
0585       if (muon.nValidRecoMuonHits() != expectedMuonIntegralValues_[4] + iOffset) {
0586         throwWithMessage("analyzeMuons, nValidRecoMuonHits does not equal expected value");
0587       }
0588       if (muon.nRecoMuonChambers() != expectedMuonIntegralValues_[5] + iOffset) {
0589         throwWithMessage("analyzeMuons, nRecoMuonChambers does not equal expected value");
0590       }
0591       if (muon.nRecoMuonChambersCSCorDT() != expectedMuonIntegralValues_[6] + iOffset) {
0592         throwWithMessage("analyzeMuons, nRecoMuonChambersCSCorDT does not equal expected value");
0593       }
0594       if (muon.nRecoMuonMatches() != expectedMuonIntegralValues_[7] + iOffset) {
0595         throwWithMessage("analyzeMuons, nRecoMuonMatches does not equal expected value");
0596       }
0597       if (muon.nRecoMuonMatchedStations() != expectedMuonIntegralValues_[8] + iOffset) {
0598         throwWithMessage("analyzeMuons, nRecoMuonMatchedStations does not equal expected value");
0599       }
0600       if (muon.nRecoMuonExpectedMatchedStations() !=
0601           static_cast<unsigned int>(expectedMuonIntegralValues_[9] + iOffset)) {
0602         throwWithMessage("analyzeMuons, nRecoMuonExpectedMatchedStations does not equal expected value");
0603       }
0604       if (muon.recoMuonStationMask() != static_cast<unsigned int>(expectedMuonIntegralValues_[10] + iOffset)) {
0605         throwWithMessage("analyzeMuons, recoMuonStationMask does not equal expected value");
0606       }
0607       if (muon.nRecoMuonMatchedRPCLayers() != expectedMuonIntegralValues_[11] + iOffset) {
0608         throwWithMessage("analyzeMuons, nRecoMuonMatchedRPCLayers does not equal expected value");
0609       }
0610       if (muon.recoMuonRPClayerMask() != static_cast<unsigned int>(expectedMuonIntegralValues_[12] + iOffset)) {
0611         throwWithMessage("analyzeMuons, recoMuonRPClayerMask does not equal expected value");
0612       }
0613       if (muon.nValidPixelHits() != expectedMuonIntegralValues_[13] + iOffset) {
0614         throwWithMessage("analyzeMuons, nValidPixelHits does not equal expected value");
0615       }
0616       if (muon.nValidStripHits() != expectedMuonIntegralValues_[14] + iOffset) {
0617         throwWithMessage("analyzeMuons, nValidStripHits does not equal expected value");
0618       }
0619       if (muon.nPixelLayersWithMeasurement() != expectedMuonIntegralValues_[15] + iOffset) {
0620         throwWithMessage("analyzeMuons, nPixelLayersWithMeasurement does not equal expected value");
0621       }
0622       if (muon.nTrackerLayersWithMeasurement() != expectedMuonIntegralValues_[16] + iOffset) {
0623         throwWithMessage("analyzeMuons, nTrackerLayersWithMeasurement does not equal expected value");
0624       }
0625       if (muon.trk_chi2() != expectedMuonFloatingPointValues_[8] + offset) {
0626         throwWithMessage("analyzeMuons, trk_chi2  does not equal expected value");
0627       }
0628       if (muon.trk_ndof() != expectedMuonFloatingPointValues_[9] + offset) {
0629         throwWithMessage("analyzeMuons, trk_ndof does not equal expected value");
0630       }
0631       if (muon.trk_dxy() != expectedMuonFloatingPointValues_[10] + offset) {
0632         throwWithMessage("analyzeMuons, trk_dxy does not equal expected value");
0633       }
0634       if (muon.trk_dz() != expectedMuonFloatingPointValues_[11] + offset) {
0635         throwWithMessage("analyzeMuons, trk_dz does not equal expected value");
0636       }
0637       if (muon.trk_qoverp() != expectedMuonFloatingPointValues_[12] + offset) {
0638         throwWithMessage("analyzeMuons, trk_qoverp does not equal expected value");
0639       }
0640       if (muon.trk_lambda() != expectedMuonFloatingPointValues_[13] + offset) {
0641         throwWithMessage("analyzeMuons, trk_lambda does not equal expected value");
0642       }
0643       if (muon.trk_pt() != expectedMuonFloatingPointValues_[14] + offset) {
0644         throwWithMessage("analyzeMuons, trk_pt does not equal expected value");
0645       }
0646       if (muon.trk_phi() != expectedMuonFloatingPointValues_[15] + offset) {
0647         throwWithMessage("analyzeMuons, trk_phi does not equal expected value");
0648       }
0649       if (muon.trk_eta() != expectedMuonFloatingPointValues_[16] + offset) {
0650         throwWithMessage("analyzeMuons, trk_eta does not equal expected value");
0651       }
0652       if (muon.trk_dxyError() != expectedMuonFloatingPointValues_[17] + offset) {
0653         throwWithMessage("analyzeMuons, trk_dxyError does not equal expected value");
0654       }
0655       if (muon.trk_dzError() != expectedMuonFloatingPointValues_[18] + offset) {
0656         throwWithMessage("analyzeMuons, trk_dzError does not equal expected value");
0657       }
0658       if (muon.trk_qoverpError() != expectedMuonFloatingPointValues_[19] + offset) {
0659         throwWithMessage("analyzeMuons, trk_qoverpError does not equal expected value");
0660       }
0661       if (muon.trk_lambdaError() != expectedMuonFloatingPointValues_[20] + offset) {
0662         throwWithMessage("analyzeMuons, trk_lambdaError does not equal expected value");
0663       }
0664       if (muon.trk_phiError() != expectedMuonFloatingPointValues_[21] + offset) {
0665         throwWithMessage("analyzeMuons, trk_phiError does not equal expected value");
0666       }
0667       if (muon.trk_dsz() != expectedMuonFloatingPointValues_[22] + offset) {
0668         throwWithMessage("analyzeMuons, trk_dsz does not equal expected value");
0669       }
0670       if (muon.trk_dszError() != expectedMuonFloatingPointValues_[23] + offset) {
0671         throwWithMessage("analyzeMuons, trk_dszError does not equal expected value");
0672       }
0673       if (muon.trk_qoverp_lambda_cov() != expectedMuonFloatingPointValues_[24] + offset) {
0674         throwWithMessage("analyzeMuons, trk_qoverp_lambda_cov does not equal expected value");
0675       }
0676       if (muon.trk_qoverp_phi_cov() != expectedMuonFloatingPointValues_[25] + offset) {
0677         throwWithMessage("analyzeMuons, trk_qoverp_phi_cov does not equal expected value");
0678       }
0679       if (muon.trk_qoverp_dxy_cov() != expectedMuonFloatingPointValues_[26] + offset) {
0680         throwWithMessage("analyzeMuons, trk_qoverp_dxy_cov does not equal expected value");
0681       }
0682       if (muon.trk_qoverp_dsz_cov() != expectedMuonFloatingPointValues_[27] + offset) {
0683         throwWithMessage("analyzeMuons, trk_qoverp_dsz_cov does not equal expected value");
0684       }
0685       if (muon.trk_lambda_phi_cov() != expectedMuonFloatingPointValues_[28] + offset) {
0686         throwWithMessage("analyzeMuons, trk_lambda_phi_cov does not equal expected value");
0687       }
0688       if (muon.trk_lambda_dxy_cov() != expectedMuonFloatingPointValues_[29] + offset) {
0689         throwWithMessage("analyzeMuons, trk_lambda_dxy_cov  does not equal expected value");
0690       }
0691       if (muon.trk_lambda_dsz_cov() != expectedMuonFloatingPointValues_[30] + offset) {
0692         throwWithMessage("analyzeMuons, trk_lambda_dsz_cov  does not equal expected value");
0693       }
0694       if (muon.trk_phi_dxy_cov() != expectedMuonFloatingPointValues_[31] + offset) {
0695         throwWithMessage("analyzeMuons, trk_phi_dxy_cov does not equal expected value");
0696       }
0697       if (muon.trk_phi_dsz_cov() != expectedMuonFloatingPointValues_[32] + offset) {
0698         throwWithMessage("analyzeMuons, trk_phi_dsz_cov does not equal expected value");
0699       }
0700       if (muon.trk_dxy_dsz_cov() != expectedMuonFloatingPointValues_[33] + offset) {
0701         throwWithMessage("analyzeMuons, trk_dxy_dsz_cov does not equal expected value");
0702       }
0703       if (muon.trk_vx() != expectedMuonFloatingPointValues_[34] + offset) {
0704         throwWithMessage("analyzeMuons, trk_vx does not equal expected value");
0705       }
0706       if (muon.trk_vy() != expectedMuonFloatingPointValues_[35] + offset) {
0707         throwWithMessage("analyzeMuons, trk_vy does not equal expected value");
0708       }
0709       if (muon.trk_vz() != expectedMuonFloatingPointValues_[36] + offset) {
0710         throwWithMessage("analyzeMuons, trk_vz does not equal expected value");
0711       }
0712       int j = 0;
0713       for (auto const& val : muon.vtxIndx()) {
0714         if (val != expectedMuonIntegralValues_[17] + iOffset + 10 * j) {
0715           throwWithMessage("analyzeMuons, vtxIndx does not contain expected value");
0716         }
0717         ++j;
0718       }
0719       if (muon.trk_hitPattern().hitCount != static_cast<uint8_t>(expectedMuonIntegralValues_[18] + iOffset)) {
0720         throwWithMessage("analyzeMuons, hitCount does not equal expected value");
0721       }
0722       if (muon.trk_hitPattern().beginTrackHits != static_cast<uint8_t>(expectedMuonIntegralValues_[19] + iOffset)) {
0723         throwWithMessage("analyzeMuons, beginTrackHits does not equal expected value");
0724       }
0725       if (muon.trk_hitPattern().endTrackHits != static_cast<uint8_t>(expectedMuonIntegralValues_[20] + iOffset)) {
0726         throwWithMessage("analyzeMuons, endTrackHits does not equal expected value");
0727       }
0728       if (muon.trk_hitPattern().beginInner != static_cast<uint8_t>(expectedMuonIntegralValues_[21] + iOffset)) {
0729         throwWithMessage("analyzeMuons, beginInner does not equal expected value");
0730       }
0731       if (muon.trk_hitPattern().endInner != static_cast<uint8_t>(expectedMuonIntegralValues_[22] + iOffset)) {
0732         throwWithMessage("analyzeMuons, endInner does not equal expected value");
0733       }
0734       if (muon.trk_hitPattern().beginOuter != static_cast<uint8_t>(expectedMuonIntegralValues_[23] + iOffset)) {
0735         throwWithMessage("analyzeMuons, beginOuter does not equal expected value");
0736       }
0737       if (muon.trk_hitPattern().endOuter != static_cast<uint8_t>(expectedMuonIntegralValues_[24] + iOffset)) {
0738         throwWithMessage("analyzeMuons, endOuter does not equal expected value");
0739       }
0740       j = 0;
0741       for (auto const& val : muon.trk_hitPattern().hitPattern) {
0742         if (val != static_cast<uint16_t>(expectedMuonIntegralValues_[25] + iOffset + 10 * j)) {
0743           throwWithMessage("analyzeMuons, hitPattern does not contain expected value");
0744         }
0745         ++j;
0746       }
0747       ++i;
0748     }
0749   }
0750 
0751   void TestReadRun3Scouting::analyzeParticles(edm::Event const& iEvent) const {
0752     auto const& particles = iEvent.get(particlesToken_);
0753     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0754     if (particles.size() != vectorSize) {
0755       throwWithMessage("analyzeParticles, particles does not have expected size");
0756     }
0757     unsigned int i = 0;
0758     for (auto const& particle : particles) {
0759       double offset = static_cast<double>(iEvent.id().event() + i);
0760       int iOffset = static_cast<int>(iEvent.id().event() + i);
0761 
0762       if (particle.pt() != expectedParticleFloatingPointValues_[0] + offset) {
0763         throwWithMessage("analyzeParticles, pt does not equal expected value");
0764       }
0765       if (particle.eta() != expectedParticleFloatingPointValues_[1] + offset) {
0766         throwWithMessage("analyzeParticles, eta does not equal expected value");
0767       }
0768       if (particle.phi() != expectedParticleFloatingPointValues_[2] + offset) {
0769         throwWithMessage("analyzeParticles, phi does not equal expected value");
0770       }
0771       if (particle.pdgId() != expectedParticleIntegralValues_[0] + iOffset) {
0772         throwWithMessage("analyzeParticles, pdgId does not equal expected value");
0773       }
0774       if (particle.vertex() != expectedParticleIntegralValues_[1] + iOffset) {
0775         throwWithMessage("analyzeParticles, vertex does not equal expected value");
0776       }
0777       if (particle.normchi2() != expectedParticleFloatingPointValues_[3] + offset) {
0778         throwWithMessage("analyzeParticles, normchi2 does not equal expected value");
0779       }
0780       if (particle.dz() != expectedParticleFloatingPointValues_[4] + offset) {
0781         throwWithMessage("analyzeParticles, dz does not equal expected value");
0782       }
0783       if (particle.dxy() != expectedParticleFloatingPointValues_[5] + offset) {
0784         throwWithMessage("analyzeParticles, dxy does not equal expected value");
0785       }
0786       if (particle.dzsig() != expectedParticleFloatingPointValues_[6] + offset) {
0787         throwWithMessage("analyzeParticles, dzsig does not equal expected value");
0788       }
0789       if (particle.dxysig() != expectedParticleFloatingPointValues_[7] + offset) {
0790         throwWithMessage("analyzeParticles, dxysig does not equal expected value");
0791       }
0792       if (particle.lostInnerHits() != static_cast<uint8_t>(expectedParticleIntegralValues_[2] + iOffset)) {
0793         throwWithMessage("analyzeParticles, lostInnerHits does not equal expected value");
0794       }
0795       if (particle.quality() != static_cast<uint8_t>(expectedParticleIntegralValues_[3] + iOffset)) {
0796         throwWithMessage("analyzeParticles, quality does not equal expected value");
0797       }
0798       if (particle.trk_pt() != expectedParticleFloatingPointValues_[8] + offset) {
0799         throwWithMessage("analyzeParticles, trk_pt does not equal expected value");
0800       }
0801       if (particle.trk_eta() != expectedParticleFloatingPointValues_[9] + offset) {
0802         throwWithMessage("analyzeParticles, trk_eta does not equal expected value");
0803       }
0804       if (particle.trk_phi() != expectedParticleFloatingPointValues_[10] + offset) {
0805         throwWithMessage("analyzeParticles, trk_phi does not equal expected value");
0806       }
0807       if (particle.relative_trk_vars() != static_cast<bool>((expectedParticleIntegralValues_[4] + iOffset) % 2)) {
0808         throwWithMessage("analyzeParticles, relative_trk_vars does not equal expected value");
0809       }
0810       ++i;
0811     }
0812   }
0813 
0814   void TestReadRun3Scouting::analyzePFJets(edm::Event const& iEvent) const {
0815     auto const& pfJets = iEvent.get(pfJetsToken_);
0816     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0817     if (pfJets.size() != vectorSize) {
0818       throwWithMessage("analyzePFJets, pfJets does not have expected size");
0819     }
0820     unsigned int i = 0;
0821     for (auto const& pfJet : pfJets) {
0822       double offset = static_cast<double>(iEvent.id().event() + i);
0823       int iOffset = static_cast<int>(iEvent.id().event() + i);
0824 
0825       if (pfJet.pt() != expectedPFJetFloatingPointValues_[0] + offset) {
0826         throwWithMessage("analyzePFJets, pt does not equal expected value");
0827       }
0828       if (pfJet.eta() != expectedPFJetFloatingPointValues_[1] + offset) {
0829         throwWithMessage("analyzePFJets, eta does not equal expected value");
0830       }
0831       if (pfJet.phi() != expectedPFJetFloatingPointValues_[2] + offset) {
0832         throwWithMessage("analyzePFJets, phi does not equal expected value");
0833       }
0834       if (pfJet.m() != expectedPFJetFloatingPointValues_[3] + offset) {
0835         throwWithMessage("analyzePFJets, m does not equal expected value");
0836       }
0837       if (pfJet.jetArea() != expectedPFJetFloatingPointValues_[4] + offset) {
0838         throwWithMessage("analyzePFJets, jetArea does not equal expected value");
0839       }
0840       if (pfJet.chargedHadronEnergy() != expectedPFJetFloatingPointValues_[5] + offset) {
0841         throwWithMessage("analyzePFJets, chargedHadronEnergy does not equal expected value");
0842       }
0843       if (pfJet.neutralHadronEnergy() != expectedPFJetFloatingPointValues_[6] + offset) {
0844         throwWithMessage("analyzePFJets, neutralHadronEnergy does not equal expected value");
0845       }
0846       if (pfJet.photonEnergy() != expectedPFJetFloatingPointValues_[7] + offset) {
0847         throwWithMessage("analyzePFJets, photonEnergy does not equal expected value");
0848       }
0849       if (pfJet.electronEnergy() != expectedPFJetFloatingPointValues_[8] + offset) {
0850         throwWithMessage("analyzePFJets, electronEnergy does not equal expected value");
0851       }
0852       if (pfJet.muonEnergy() != expectedPFJetFloatingPointValues_[9] + offset) {
0853         throwWithMessage("analyzePFJets, muonEnergy does not equal expected value");
0854       }
0855       if (pfJet.HFHadronEnergy() != expectedPFJetFloatingPointValues_[10] + offset) {
0856         throwWithMessage("analyzePFJets, HFHadronEnergy does not equal expected value");
0857       }
0858       if (pfJet.HFEMEnergy() != expectedPFJetFloatingPointValues_[11] + offset) {
0859         throwWithMessage("analyzePFJets, HFEMEnergy does not equal expected value");
0860       }
0861       if (pfJet.chargedHadronMultiplicity() != expectedPFJetIntegralValues_[0] + iOffset) {
0862         throwWithMessage("analyzePFJets, chargedHadronMultiplicity does not equal expected value");
0863       }
0864       if (pfJet.neutralHadronMultiplicity() != expectedPFJetIntegralValues_[1] + iOffset) {
0865         throwWithMessage("analyzePFJets, neutralHadronMultiplicity does not equal expected value");
0866       }
0867       if (pfJet.photonMultiplicity() != expectedPFJetIntegralValues_[2] + iOffset) {
0868         throwWithMessage("analyzePFJets, photonMultiplicity does not equal expected value");
0869       }
0870       if (pfJet.electronMultiplicity() != expectedPFJetIntegralValues_[3] + iOffset) {
0871         throwWithMessage("analyzePFJets, electronMultiplicity does not equal expected value");
0872       }
0873       if (pfJet.muonMultiplicity() != expectedPFJetIntegralValues_[4] + iOffset) {
0874         throwWithMessage("analyzePFJets, muonMultiplicity does not equal expected value");
0875       }
0876       if (pfJet.HFHadronMultiplicity() != expectedPFJetIntegralValues_[5] + iOffset) {
0877         throwWithMessage("analyzePFJets, HFHadronMultiplicity does not equal expected value");
0878       }
0879       if (pfJet.HFEMMultiplicity() != expectedPFJetIntegralValues_[6] + iOffset) {
0880         throwWithMessage("analyzePFJets, HFEMMultiplicity does not equal expected value");
0881       }
0882       if (pfJet.HOEnergy() != expectedPFJetFloatingPointValues_[12] + offset) {
0883         throwWithMessage("analyzePFJets, HOEnergy does not equal expected value");
0884       }
0885       if (pfJet.csv() != expectedPFJetFloatingPointValues_[13] + offset) {
0886         throwWithMessage("analyzePFJets, csv does not equal expected value");
0887       }
0888       if (pfJet.mvaDiscriminator() != expectedPFJetFloatingPointValues_[14] + offset) {
0889         throwWithMessage("analyzePFJets, mvaDiscriminator does not equal expected value");
0890       }
0891       int j = 0;
0892       for (auto const& val : pfJet.constituents()) {
0893         if (val != expectedPFJetIntegralValues_[7] + iOffset + 10 * j) {
0894           throwWithMessage("analyzePFJets, constituents does not contain expected value");
0895         }
0896         ++j;
0897       }
0898       ++i;
0899     }
0900   }
0901 
0902   void TestReadRun3Scouting::analyzePhotons(edm::Event const& iEvent) const {
0903     auto const& photons = iEvent.get(photonsToken_);
0904     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0905     if (photons.size() != vectorSize) {
0906       throwWithMessage("analyzePhotons, photons does not have expected size");
0907     }
0908     unsigned int i = 0;
0909     for (auto const& photon : photons) {
0910       double offset = static_cast<double>(iEvent.id().event() + i);
0911       int iOffset = static_cast<int>(iEvent.id().event() + i);
0912 
0913       if (photon.pt() != expectedPhotonFloatingPointValues_[0] + offset) {
0914         throwWithMessage("analyzePhotons, pt does not equal expected value");
0915       }
0916       if (photon.eta() != expectedPhotonFloatingPointValues_[1] + offset) {
0917         throwWithMessage("analyzePhotons, eta does not equal expected value");
0918       }
0919       if (photon.phi() != expectedPhotonFloatingPointValues_[2] + offset) {
0920         throwWithMessage("analyzePhotons, phi does not equal expected value");
0921       }
0922       if (photon.m() != expectedPhotonFloatingPointValues_[3] + offset) {
0923         throwWithMessage("analyzePhotons, m does not equal expected value");
0924       }
0925       if (photon.sigmaIetaIeta() != expectedPhotonFloatingPointValues_[4] + offset) {
0926         throwWithMessage("analyzePhotons, sigmaIetaIeta does not equal expected value");
0927       }
0928       if (photon.hOverE() != expectedPhotonFloatingPointValues_[5] + offset) {
0929         throwWithMessage("analyzePhotons, hOverE does not equal expected value");
0930       }
0931       if (photon.ecalIso() != expectedPhotonFloatingPointValues_[6] + offset) {
0932         throwWithMessage("analyzePhotons, ecalIso does not equal expected value");
0933       }
0934       if (photon.hcalIso() != expectedPhotonFloatingPointValues_[7] + offset) {
0935         throwWithMessage("analyzePhotons, hcalIso does not equal expected value");
0936       }
0937       if (photon.trkIso() != expectedPhotonFloatingPointValues_[8] + offset) {
0938         throwWithMessage("analyzePhotons, trkIso does not equal expected value");
0939       }
0940       if (photon.r9() != expectedPhotonFloatingPointValues_[9] + offset) {
0941         throwWithMessage("analyzePhotons, r9 does not equal expected value");
0942       }
0943       if (photon.sMin() != expectedPhotonFloatingPointValues_[10] + offset) {
0944         throwWithMessage("analyzePhotons, sMin does not equal expected value");
0945       }
0946       if (photon.sMaj() != expectedPhotonFloatingPointValues_[11] + offset) {
0947         throwWithMessage("analyzePhotons, sMaj does not equal expected value");
0948       }
0949       if (photon.seedId() != static_cast<unsigned int>(expectedPhotonIntegralValues_[0] + iOffset)) {
0950         throwWithMessage("analyzePhotons, seedId does not equal expected value");
0951       }
0952 
0953       if (photon.energyMatrix().size() != vectorSize) {
0954         throwWithMessage("analyzePhotons, energyMatrix does not have expected size");
0955       }
0956       unsigned int j = 0;
0957       for (auto const& val : photon.energyMatrix()) {
0958         if (val != expectedPhotonFloatingPointValues_[12] + offset + 10 * j) {
0959           throwWithMessage("analyzePhotons, energyMatrix does not contain expected value");
0960         }
0961         ++j;
0962       }
0963       if (photon.detIds().size() != vectorSize) {
0964         throwWithMessage("analyzePhotons, detIds does not have expected size");
0965       }
0966       j = 0;
0967       for (auto const& val : photon.detIds()) {
0968         if (val != static_cast<uint32_t>(expectedPhotonIntegralValues_[1] + iOffset + 10 * j)) {
0969           throwWithMessage("analyzePhotons, detIds does not contain expected value");
0970         }
0971         ++j;
0972       }
0973       if (photon.timingMatrix().size() != vectorSize) {
0974         throwWithMessage("analyzePhotons, timingMatrix does not have expected size");
0975       }
0976       j = 0;
0977       for (auto const& val : photon.timingMatrix()) {
0978         if (val != expectedPhotonFloatingPointValues_[13] + offset + 10 * j) {
0979           throwWithMessage("analyzePhotons, timingMatrix does not contain expected value");
0980         }
0981         ++j;
0982       }
0983       if (photon.rechitZeroSuppression() != static_cast<bool>((expectedPhotonIntegralValues_[2] + iOffset) % 2)) {
0984         throwWithMessage("analyzePhotons, rechitZeroSuppression does not equal expected value");
0985       }
0986       if (inputPhotonClassVersion_ == 6) {
0987         if (photon.rawEnergy() != expectedPhotonFloatingPointValues_[14] + offset) {
0988           throwWithMessage("analyzePhotons, rawEnergy does not equal expected value");
0989         }
0990         if (photon.preshowerEnergy() != expectedPhotonFloatingPointValues_[15] + offset) {
0991           throwWithMessage("analyzePhotons, preshowerEnergy does not equal expected value");
0992         }
0993         if (photon.corrEcalEnergyError() != expectedPhotonFloatingPointValues_[16] + offset) {
0994           throwWithMessage("analyzePhotons, corrEcalEnergyError does not equal expected value");
0995         }
0996         if (photon.nClusters() != static_cast<unsigned int>(expectedPhotonIntegralValues_[3] + iOffset)) {
0997           throwWithMessage("analyzePhotons, nClusters does not equal expected value");
0998         }
0999         if (photon.nCrystals() != static_cast<unsigned int>(expectedPhotonIntegralValues_[4] + iOffset)) {
1000           throwWithMessage("analyzePhotons, nCrystals does not equal expected value");
1001         }
1002       }
1003       ++i;
1004     }
1005   }
1006 
1007   void TestReadRun3Scouting::analyzeTracks(edm::Event const& iEvent) const {
1008     auto const& tracks = iEvent.get(tracksToken_);
1009     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
1010     if (tracks.size() != vectorSize) {
1011       throwWithMessage("analyzeTracks, tracks does not have expected size");
1012     }
1013     unsigned int i = 0;
1014     for (auto const& track : tracks) {
1015       double offset = static_cast<double>(iEvent.id().event() + i);
1016       int iOffset = static_cast<int>(iEvent.id().event() + i);
1017 
1018       if (track.tk_pt() != expectedTrackFloatingPointValues_[0] + offset) {
1019         throwWithMessage("analyzeTracks, tk_pt does not equal expected value");
1020       }
1021       if (track.tk_eta() != expectedTrackFloatingPointValues_[1] + offset) {
1022         throwWithMessage("analyzeTracks, tk_eta does not equal expected value");
1023       }
1024       if (track.tk_phi() != expectedTrackFloatingPointValues_[2] + offset) {
1025         throwWithMessage("analyzeTracks, tk_phi does not equal expected value");
1026       }
1027       if (track.tk_chi2() != expectedTrackFloatingPointValues_[3] + offset) {
1028         throwWithMessage("analyzeTracks, tk_chi2 does not equal expected value");
1029       }
1030       if (track.tk_ndof() != expectedTrackFloatingPointValues_[4] + offset) {
1031         throwWithMessage("analyzeTracks, tk_ndof does not equal expected value");
1032       }
1033       if (track.tk_charge() != expectedTrackIntegralValues_[0] + iOffset) {
1034         throwWithMessage("analyzeTracks, tk_charge does not equal expected value");
1035       }
1036       if (track.tk_dxy() != expectedTrackFloatingPointValues_[5] + offset) {
1037         throwWithMessage("analyzeTracks, tk_dxy does not equal expected value");
1038       }
1039       if (track.tk_dz() != expectedTrackFloatingPointValues_[6] + offset) {
1040         throwWithMessage("analyzeTracks, tk_dz does not equal expected value");
1041       }
1042       if (track.tk_nValidPixelHits() != expectedTrackIntegralValues_[1] + iOffset) {
1043         throwWithMessage("analyzeTracks, tk_nValidPixelHits does not equal expected value");
1044       }
1045       if (track.tk_nTrackerLayersWithMeasurement() != expectedTrackIntegralValues_[2] + iOffset) {
1046         throwWithMessage("analyzeTracks, tk_nTrackerLayersWithMeasurement does not equal expected value");
1047       }
1048       if (track.tk_nValidStripHits() != expectedTrackIntegralValues_[3] + iOffset) {
1049         throwWithMessage("analyzeTracks, tk_nValidStripHits does not equal expected value");
1050       }
1051       if (track.tk_qoverp() != expectedTrackFloatingPointValues_[7] + offset) {
1052         throwWithMessage("analyzeTracks, tk_qoverp does not equal expected value");
1053       }
1054       if (track.tk_lambda() != expectedTrackFloatingPointValues_[8] + offset) {
1055         throwWithMessage("analyzeTracks, tk_lambda does not equal expected value");
1056       }
1057       if (track.tk_dxy_Error() != expectedTrackFloatingPointValues_[9] + offset) {
1058         throwWithMessage("analyzeTracks, tk_dxy_Error does not equal expected value");
1059       }
1060       if (track.tk_dz_Error() != expectedTrackFloatingPointValues_[10] + offset) {
1061         throwWithMessage("analyzeTracks, tk_dz_Error does not equal expected value");
1062       }
1063       if (track.tk_qoverp_Error() != expectedTrackFloatingPointValues_[11] + offset) {
1064         throwWithMessage("analyzeTracks, tk_qoverp_Error does not equal expected value");
1065       }
1066       if (track.tk_lambda_Error() != expectedTrackFloatingPointValues_[12] + offset) {
1067         throwWithMessage("analyzeTracks, tk_lambda_Error does not equal expected value");
1068       }
1069       if (track.tk_phi_Error() != expectedTrackFloatingPointValues_[13] + offset) {
1070         throwWithMessage("analyzeTracks, tk_phi_Error does not equal expected value");
1071       }
1072       if (track.tk_dsz() != expectedTrackFloatingPointValues_[14] + offset) {
1073         throwWithMessage("analyzeTracks, tk_dsz does not equal expected value");
1074       }
1075       if (track.tk_dsz_Error() != expectedTrackFloatingPointValues_[15] + offset) {
1076         throwWithMessage("analyzeTracks, tk_dsz_Error does not equal expected value");
1077       }
1078       if (track.tk_qoverp_lambda_cov() != expectedTrackFloatingPointValues_[16] + offset) {
1079         throwWithMessage("analyzeTracks, tk_qoverp_lambda_cov does not equal expected value");
1080       }
1081       if (track.tk_qoverp_phi_cov() != expectedTrackFloatingPointValues_[17] + offset) {
1082         throwWithMessage("analyzeTracks, tk_qoverp_phi_cov does not equal expected value");
1083       }
1084       if (track.tk_qoverp_dxy_cov() != expectedTrackFloatingPointValues_[18] + offset) {
1085         throwWithMessage("analyzeTracks, tk_qoverp_dxy_cov does not equal expected value");
1086       }
1087       if (track.tk_qoverp_dsz_cov() != expectedTrackFloatingPointValues_[19] + offset) {
1088         throwWithMessage("analyzeTracks, tk_qoverp_dsz_cov does not equal expected value");
1089       }
1090       if (track.tk_lambda_phi_cov() != expectedTrackFloatingPointValues_[20] + offset) {
1091         throwWithMessage("analyzeTracks, tk_lambda_phi_cov does not equal expected value");
1092       }
1093       if (track.tk_lambda_dxy_cov() != expectedTrackFloatingPointValues_[21] + offset) {
1094         throwWithMessage("analyzeTracks, tk_lambda_dxy_cov does not equal expected value");
1095       }
1096       if (track.tk_lambda_dsz_cov() != expectedTrackFloatingPointValues_[22] + offset) {
1097         throwWithMessage("analyzeTracks, tk_lambda_dsz_cov does not equal expected value");
1098       }
1099       if (track.tk_phi_dxy_cov() != expectedTrackFloatingPointValues_[23] + offset) {
1100         throwWithMessage("analyzeTracks, tk_phi_dxy_cov does not equal expected value");
1101       }
1102       if (track.tk_phi_dsz_cov() != expectedTrackFloatingPointValues_[24] + offset) {
1103         throwWithMessage("analyzeTracks, tk_phi_dsz_cov does not equal expected value");
1104       }
1105       if (track.tk_dxy_dsz_cov() != expectedTrackFloatingPointValues_[25] + offset) {
1106         throwWithMessage("analyzeTracks, tk_dxy_dsz_cov does not equal expected value");
1107       }
1108       if (track.tk_vtxInd() != expectedTrackIntegralValues_[4] + iOffset) {
1109         throwWithMessage("analyzeTracks, tk_vtxInd does not equal expected value");
1110       }
1111       if (track.tk_vx() != expectedTrackFloatingPointValues_[26] + offset) {
1112         throwWithMessage("analyzeTracks, tk_vx does not equal expected value");
1113       }
1114       if (track.tk_vy() != expectedTrackFloatingPointValues_[27] + offset) {
1115         throwWithMessage("analyzeTracks, tk_vy does not equal expected value");
1116       }
1117       if (track.tk_vz() != expectedTrackFloatingPointValues_[28] + offset) {
1118         throwWithMessage("analyzeTracks, float tk_vz does not equal expected value");
1119       }
1120       ++i;
1121     }
1122   }
1123 
1124   void TestReadRun3Scouting::analyzeVertexes(edm::Event const& iEvent) const {
1125     auto const& vertexes = iEvent.get(vertexesToken_);
1126     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
1127     if (vertexes.size() != vectorSize) {
1128       throwWithMessage("analyzeVertexes, vertexes does not have expected size");
1129     }
1130     unsigned int i = 0;
1131     for (auto const& vertex : vertexes) {
1132       double offset = static_cast<double>(iEvent.id().event() + i);
1133       int iOffset = static_cast<int>(iEvent.id().event() + i);
1134 
1135       if (vertex.x() != expectedVertexFloatingPointValues_[0] + offset) {
1136         throwWithMessage("analyzeVertexes, x does not equal expected value");
1137       }
1138       if (vertex.y() != expectedVertexFloatingPointValues_[1] + offset) {
1139         throwWithMessage("analyzeVertexes, y does not equal expected value");
1140       }
1141       if (vertex.z() != expectedVertexFloatingPointValues_[2] + offset) {
1142         throwWithMessage("analyzeVertexes, z does not equal expected value");
1143       }
1144       if (vertex.zError() != expectedVertexFloatingPointValues_[3] + offset) {
1145         throwWithMessage("analyzeVertexes, zError does not equal expected value");
1146       }
1147       if (vertex.xError() != expectedVertexFloatingPointValues_[4] + offset) {
1148         throwWithMessage("analyzeVertexes, xError does not equal expected value");
1149       }
1150       if (vertex.yError() != expectedVertexFloatingPointValues_[5] + offset) {
1151         throwWithMessage("analyzeVertexes, yError does not equal expected value");
1152       }
1153       if (vertex.tracksSize() != expectedVertexIntegralValues_[0] + iOffset) {
1154         throwWithMessage("analyzeVertexes, tracksSize does not equal expected value");
1155       }
1156       if (vertex.chi2() != expectedVertexFloatingPointValues_[6] + offset) {
1157         throwWithMessage("analyzeVertexes, chi2 does not equal expected value");
1158       }
1159       if (vertex.ndof() != expectedVertexIntegralValues_[1] + iOffset) {
1160         throwWithMessage("analyzeVertexes, ndof does not equal expected value");
1161       }
1162       if (vertex.isValidVtx() != static_cast<bool>((expectedVertexIntegralValues_[2] + iOffset) % 2)) {
1163         throwWithMessage("analyzeVertexes, isValidVtx does not equal expected value");
1164       }
1165 
1166       if (inputVertexClassVersion_ == 4) {
1167         if (vertex.xyCov() != expectedVertexFloatingPointValues_[7] + offset) {
1168           throwWithMessage("analyzeVertexes, xy cov. does not equal expected value");
1169         }
1170         if (vertex.xzCov() != expectedVertexFloatingPointValues_[8] + offset) {
1171           throwWithMessage("analyzeVertexes, xz cov. does not equal expected value");
1172         }
1173         if (vertex.yzCov() != expectedVertexFloatingPointValues_[9] + offset) {
1174           throwWithMessage("analyzeVertexes, yz cov. does not equal expected value");
1175         }
1176       }
1177       ++i;
1178     }
1179   }
1180 
1181   void TestReadRun3Scouting::throwWithMessageFromConstructor(const char* msg) const {
1182     throw cms::Exception("TestFailure") << "TestReadRun3Scouting constructor, " << msg;
1183   }
1184 
1185   void TestReadRun3Scouting::throwWithMessage(const char* msg) const {
1186     throw cms::Exception("TestFailure") << "TestReadRun3Scouting::analyze, " << msg;
1187   }
1188 
1189 }  // namespace edmtest
1190 
1191 using edmtest::TestReadRun3Scouting;
1192 DEFINE_FWK_MODULE(TestReadRun3Scouting);