Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:08

0001 // -*- C++ -*-
0002 //
0003 // Package:    IOPool/Input
0004 // Class:      SchemaEvolutionTestRead
0005 //
0006 /**\class edmtest::SchemaEvolutionTestRead
0007   Description: Used as part of tests of ROOT's schema evolution
0008   features. These features allow reading a persistent object
0009   whose data format has changed since it was written.
0010 */
0011 // Original Author:  W. David Dagenhart
0012 //         Created:  28 July 2023
0013 
0014 #include "DataFormats/TestObjects/interface/VectorVectorTop.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/Frameworkfwd.h"
0017 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0022 #include "FWCore/Utilities/interface/EDGetToken.h"
0023 #include "FWCore/Utilities/interface/Exception.h"
0024 #include "FWCore/Utilities/interface/InputTag.h"
0025 #include "FWCore/Utilities/interface/StreamID.h"
0026 
0027 #include <vector>
0028 
0029 namespace edmtest {
0030 
0031   class SchemaEvolutionTestRead : public edm::global::EDAnalyzer<> {
0032   public:
0033     SchemaEvolutionTestRead(edm::ParameterSet const&);
0034     void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;
0035     static void fillDescriptions(edm::ConfigurationDescriptions&);
0036 
0037   private:
0038     void analyzeVectorVector(edm::Event const&) const;
0039     void analyzeVectorVectorNonSplit(edm::Event const&) const;
0040 
0041     void throwWithMessageFromConstructor(const char*) const;
0042     void throwWithMessage(const char*) const;
0043 
0044     // These expected values are meaningless other than we use them
0045     // to check that values read from persistent storage match the values
0046     // we know were written.
0047 
0048     const std::vector<int> expectedVectorVectorIntegralValues_;
0049     const edm::EDGetTokenT<VectorVectorTop> vectorVectorToken_;
0050     const edm::EDGetTokenT<VectorVectorTopNonSplit> vectorVectorNonSplitToken_;
0051   };
0052 
0053   SchemaEvolutionTestRead::SchemaEvolutionTestRead(edm::ParameterSet const& iPSet)
0054       : expectedVectorVectorIntegralValues_(iPSet.getParameter<std::vector<int>>("expectedVectorVectorIntegralValues")),
0055         vectorVectorToken_(consumes(iPSet.getParameter<edm::InputTag>("vectorVectorTag"))),
0056         vectorVectorNonSplitToken_(consumes(iPSet.getParameter<edm::InputTag>("vectorVectorTag"))) {
0057     if (expectedVectorVectorIntegralValues_.size() != 15) {
0058       throwWithMessageFromConstructor("test configuration error, expectedVectorVectorIntegralValues must have size 15");
0059     }
0060   }
0061 
0062   void SchemaEvolutionTestRead::analyze(edm::StreamID, edm::Event const& iEvent, edm::EventSetup const&) const {
0063     analyzeVectorVector(iEvent);
0064     analyzeVectorVectorNonSplit(iEvent);
0065   }
0066 
0067   void SchemaEvolutionTestRead::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0068     edm::ParameterSetDescription desc;
0069     desc.add<std::vector<int>>("expectedVectorVectorIntegralValues");
0070     desc.add<edm::InputTag>("vectorVectorTag");
0071     descriptions.addDefault(desc);
0072   }
0073 
0074   void SchemaEvolutionTestRead::analyzeVectorVector(edm::Event const& iEvent) const {
0075     auto const& vectorVector = iEvent.get(vectorVectorToken_);
0076     unsigned int vectorSize = 2 + iEvent.id().event() % 4;
0077     if (vectorVector.outerVector_.size() != vectorSize) {
0078       throwWithMessage("analyzeVectorVector, vectorVector does not have expected size");
0079     }
0080     unsigned int i = 0;
0081     for (auto const& middleVector : vectorVector.outerVector_) {
0082       if (middleVector.middleVector_.size() != vectorSize) {
0083         throwWithMessage("analyzeVectorVector, middleVector does not have expected size");
0084       }
0085       int iOffset = static_cast<int>(iEvent.id().event() + i);
0086 
0087       int j = 0;
0088       for (auto const& element : middleVector.middleVector_) {
0089         if (element.a_ != expectedVectorVectorIntegralValues_[0] + iOffset + j * 10) {
0090           throwWithMessage("analyzeVectorVector, element a_ does not contain expected value");
0091         }
0092         if (element.b_ != expectedVectorVectorIntegralValues_[0] + iOffset + j * 100) {
0093           throwWithMessage("analyzeVectorVector, element b_ does not contain expected value");
0094         }
0095 
0096         SchemaEvolutionChangeOrder const& changeOrder = element.changeOrder_;
0097         if (changeOrder.a_ != expectedVectorVectorIntegralValues_[1] + iOffset + j * 11) {
0098           throwWithMessage("analyzeVectorVector, changeOrder a_ does not contain expected value");
0099         }
0100         if (changeOrder.b_ != expectedVectorVectorIntegralValues_[1] + iOffset + j * 101) {
0101           throwWithMessage("analyzeVectorVector, changeOrder b_ does not contain expected value");
0102         }
0103 
0104         SchemaEvolutionAddMember const& addMember = element.addMember_;
0105         if (addMember.a_ != expectedVectorVectorIntegralValues_[2] + iOffset + j * 12) {
0106           throwWithMessage("analyzeVectorVector, addMember a_ does not contain expected value");
0107         }
0108         if (addMember.b_ != expectedVectorVectorIntegralValues_[2] + iOffset + j * 102) {
0109           throwWithMessage("analyzeVectorVector, addMember b_ does not contain expected value");
0110         }
0111 
0112         SchemaEvolutionRemoveMember const& removeMember = element.removeMember_;
0113         if (removeMember.a_ != expectedVectorVectorIntegralValues_[3] + iOffset + j * 13) {
0114           throwWithMessage("analyzeVectorVector, removeMember a_ does not contain expected value");
0115         }
0116 
0117         SchemaEvolutionMoveToBase const& moveToBase = element.moveToBase_;
0118         if (moveToBase.a_ != expectedVectorVectorIntegralValues_[4] + iOffset + j * 14) {
0119           throwWithMessage("analyzeVectorVector, moveToBase a_ does not contain expected value");
0120         }
0121         if (moveToBase.b_ != expectedVectorVectorIntegralValues_[4] + iOffset + j * 104) {
0122           throwWithMessage("analyzeVectorVector, moveToBase b_ does not contain expected value");
0123         }
0124         if (moveToBase.c_ != expectedVectorVectorIntegralValues_[4] + iOffset + j * 1004) {
0125           throwWithMessage("analyzeVectorVector, moveToBase c_ does not contain expected value");
0126         }
0127         if (moveToBase.d_ != expectedVectorVectorIntegralValues_[4] + iOffset + j * 10004) {
0128           throwWithMessage("analyzeVectorVector, moveToBase d_ does not contain expected value");
0129         }
0130 
0131         SchemaEvolutionChangeType const& changeType = element.changeType_;
0132         if (static_cast<int>(changeType.a_) != expectedVectorVectorIntegralValues_[5] + iOffset + j * 15) {
0133           throwWithMessage("analyzeVectorVector, changeType a_ does not contain expected value");
0134         }
0135         if (static_cast<int>(changeType.b_) != expectedVectorVectorIntegralValues_[5] + iOffset + j * 105) {
0136           throwWithMessage("analyzeVectorVector, changeType b_ does not contain expected value");
0137         }
0138 
0139         SchemaEvolutionAddBase const& addBase = element.addBase_;
0140         if (addBase.a_ != expectedVectorVectorIntegralValues_[6] + iOffset + j * 16) {
0141           throwWithMessage("analyzeVectorVector, addToBase a_ does not contain expected value");
0142         }
0143         if (addBase.b_ != expectedVectorVectorIntegralValues_[6] + iOffset + j * 106) {
0144           throwWithMessage("analyzeVectorVector, addToBase b_ does not contain expected value");
0145         }
0146 
0147         SchemaEvolutionPointerToMember const& pointerToMember = element.pointerToMember_;
0148         if (pointerToMember.a_ != expectedVectorVectorIntegralValues_[7] + iOffset + j * 17) {
0149           throwWithMessage("analyzeVectorVector, pointerToMember a_ does not contain expected value");
0150         }
0151         if (pointerToMember.b_ != expectedVectorVectorIntegralValues_[7] + iOffset + j * 107) {
0152           throwWithMessage("analyzeVectorVector, pointerToMember b_ does not contain expected value");
0153         }
0154         // This part is commented out because it fails. My conclusion is that ROOT
0155         // does not properly support schema evolution in this case. CMS does not
0156         // usually use pointers in persistent formats. So for now we are just ignoring
0157         // this issue.
0158         // if (pointerToMember.c() != expectedVectorVectorIntegralValues_[7] + iOffset + j * 1007) {
0159         //   throwWithMessage("analyzeVectorVector, pointerToMember c_ does not contain expected value");
0160         // }
0161 
0162         SchemaEvolutionPointerToUniquePtr const& pointerToUniquePtr = element.pointerToUniquePtr_;
0163         if (pointerToUniquePtr.a_ != expectedVectorVectorIntegralValues_[8] + iOffset + j * 18) {
0164           throwWithMessage("analyzeVectorVector, pointerToUniquePtr a_ does not contain expected value");
0165         }
0166         if (pointerToUniquePtr.b_ != expectedVectorVectorIntegralValues_[8] + iOffset + j * 108) {
0167           throwWithMessage("analyzeVectorVector, pointerToUniquePtr b_ does not contain expected value");
0168         }
0169         if (pointerToUniquePtr.contained_->c_ != expectedVectorVectorIntegralValues_[8] + iOffset + j * 1008) {
0170           throwWithMessage("analyzeVectorVector, pointerToUniquePtr c_ does not contain expected value");
0171         }
0172 
0173         SchemaEvolutionCArrayToStdArray const& cArrayToStdArray = element.cArrayToStdArray_;
0174         if (cArrayToStdArray.a_[0] != expectedVectorVectorIntegralValues_[9] + iOffset + j * 19) {
0175           throwWithMessage("analyzeVectorVector, cArrayToStdArray a_[0] does not contain expected value");
0176         }
0177         if (cArrayToStdArray.a_[1] != expectedVectorVectorIntegralValues_[9] + iOffset + j * 109) {
0178           throwWithMessage("analyzeVectorVector, cArrayToStdArray a_[1] does not contain expected value");
0179         }
0180         if (cArrayToStdArray.a_[2] != expectedVectorVectorIntegralValues_[9] + iOffset + j * 1009) {
0181           throwWithMessage("analyzeVectorVector, cArrayToStdArray a_[2] does not contain expected value");
0182         }
0183 
0184         // This part is commented out because it fails. My conclusion is that ROOT
0185         // does not properly support schema evolution in this case. CMS does not
0186         // usually use pointers in persistent formats. So for now we are just ignoring
0187         // this issue. Note that pointerToMember fails because the values are incorrect.
0188         // This one is also commented out of the format, the write function and here
0189         // because simply reading the object causes a fatal exception even without
0190         // checking the values.
0191         // SchemaEvolutionCArrayToStdVector const& cArrayToStdVector = element.cArrayToStdVector_;
0192         // if (cArrayToStdVector.a_[0] != expectedVectorVectorIntegralValues_[10] + iOffset + j * 20) {
0193         //   throwWithMessage("analyzeVectorVector, cArrayToStdVector a_[0] does not contain expected value");
0194         // }
0195         // if (cArrayToStdVector.a_[1] != expectedVectorVectorIntegralValues_[10] + iOffset + j * 110) {
0196         //   throwWithMessage("analyzeVectorVector, cArrayToStdVector a_[1] does not contain expected value");
0197         // }
0198         // if (cArrayToStdVector.a_[2] != expectedVectorVectorIntegralValues_[10] + iOffset + j * 1010) {
0199         //  throwWithMessage("analyzeVectorVector, cArrayToStdVector a_[2] does not contain expected value");
0200         // }
0201 
0202         {
0203           SchemaEvolutionVectorToList const& vectorToList = element.vectorToList_;
0204           auto iter = vectorToList.a_.cbegin();
0205           auto iter0 = iter;
0206           auto iter1 = ++iter;
0207           auto iter2 = ++iter;
0208           if (*iter0 != expectedVectorVectorIntegralValues_[11] + iOffset + j * 21) {
0209             throwWithMessage("vectorToList, element 0 does not contain expected value");
0210           }
0211           if (*iter1 != expectedVectorVectorIntegralValues_[11] + iOffset + j * 111) {
0212             throwWithMessage("vectorToList, element 1 does not contain expected value");
0213           }
0214           if (*iter2 != expectedVectorVectorIntegralValues_[11] + iOffset + j * 1011) {
0215             throwWithMessage("vectorToList, element 2 does not contain expected value");
0216           }
0217         }
0218 
0219         {
0220           SchemaEvolutionMapToUnorderedMap const& mapToUnorderedMap = element.mapToUnorderedMap_;
0221           if (mapToUnorderedMap.a_.size() != 3) {
0222             throwWithMessage("mapToUnorderedMap, map has unexpected size");
0223           }
0224           auto iter = mapToUnorderedMap.a_.cbegin();
0225 
0226           // Easier to check values if we sort them first so sort them in a regular map
0227           std::map<int, int> orderedMap;
0228           orderedMap.insert(*iter);
0229           ++iter;
0230           orderedMap.insert(*iter);
0231           ++iter;
0232           orderedMap.insert(*iter);
0233 
0234           auto orderedIter = orderedMap.cbegin();
0235           auto iter0 = orderedIter;
0236           auto iter1 = ++orderedIter;
0237           auto iter2 = ++orderedIter;
0238           if (iter0->first != expectedVectorVectorIntegralValues_[12] + iOffset + j * 22) {
0239             throwWithMessage("mapToUnorderedMap, element 0 key does not contain expected value");
0240           }
0241           if (iter0->second != expectedVectorVectorIntegralValues_[12] + iOffset + j * 112) {
0242             throwWithMessage("mapToUnorderedMap, element 0 does not contain expected value");
0243           }
0244           if (iter1->first != expectedVectorVectorIntegralValues_[12] + iOffset + j * 1012 + 1012) {
0245             throwWithMessage("mapToUnorderedMap, element 1 key does not contain expected value");
0246           }
0247           if (iter1->second != expectedVectorVectorIntegralValues_[12] + iOffset + j * 10012) {
0248             throwWithMessage("mapToUnorderedMap, element 1 does not contain expected value");
0249           }
0250           if (iter2->first != expectedVectorVectorIntegralValues_[12] + iOffset + j * 100012 + 100012) {
0251             throwWithMessage("mapToUnorderedMap, element 2 key does not contain expected value");
0252           }
0253           if (iter2->second != expectedVectorVectorIntegralValues_[12] + iOffset + j * 1000012) {
0254             throwWithMessage("mapToUnorderedMap, element 2 does not contain expected value");
0255           }
0256         }
0257         ++j;
0258       }
0259       ++i;
0260     }
0261   }
0262 
0263   void SchemaEvolutionTestRead::analyzeVectorVectorNonSplit(edm::Event const& iEvent) const {
0264     auto const& vectorVectorNonSplit = iEvent.get(vectorVectorNonSplitToken_);
0265     unsigned int vectorSize = 1;
0266     if (vectorVectorNonSplit.outerVector_.size() != vectorSize) {
0267       throwWithMessage("analyzeVectorVectorNonSplit, outerVector does not have expected size");
0268     }
0269     for (auto const& middleVector : vectorVectorNonSplit.outerVector_) {
0270       if (middleVector.middleVector_.size() != vectorSize) {
0271         throwWithMessage("analyzeVectorVectorNonSplit, middleVector does not have expected size");
0272       }
0273       for (auto const& element : middleVector.middleVector_) {
0274         if (element.a_ != expectedVectorVectorIntegralValues_[13]) {
0275           throwWithMessage("analyzeVectorVectorNonSplit, element a_ does not contain expected value");
0276         }
0277       }
0278     }
0279   }
0280 
0281   void SchemaEvolutionTestRead::throwWithMessageFromConstructor(const char* msg) const {
0282     throw cms::Exception("TestFailure") << "SchemaEvolutionTestRead constructor, " << msg;
0283   }
0284 
0285   void SchemaEvolutionTestRead::throwWithMessage(const char* msg) const {
0286     throw cms::Exception("TestFailure") << "SchemaEvolutionTestRead::analyze, " << msg;
0287   }
0288 
0289 }  // namespace edmtest
0290 
0291 using edmtest::SchemaEvolutionTestRead;
0292 DEFINE_FWK_MODULE(SchemaEvolutionTestRead);