Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-08 22:25:10

0001 #include <algorithm>
0002 #include <cassert>
0003 #include <vector>
0004 
0005 #include "FWCore/Sources/interface/DaqProvenanceHelper.h"
0006 
0007 #include "DataFormats/Provenance/interface/ProductDependencies.h"
0008 #include "DataFormats/Provenance/interface/BranchIDList.h"
0009 #include "DataFormats/Provenance/interface/ProcessHistory.h"
0010 #include "DataFormats/Provenance/interface/ProcessHistoryRegistry.h"
0011 #include "DataFormats/Provenance/interface/ProductRegistry.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/Utilities/interface/EDMException.h"
0014 #include "FWCore/Reflection/interface/TypeWithDict.h"
0015 #include "FWCore/Version/interface/GetReleaseVersion.h"
0016 
0017 namespace {
0018   edm::ProductDescription makeDescriptionForDaqProvHelper(edm::TypeID const& rawDataType,
0019                                                           std::string const& collectionName,
0020                                                           std::string const& friendlyName,
0021                                                           std::string const& sourceLabel) {
0022     edm::ProductDescription desc(edm::InEvent,
0023                                  "rawDataCollector",
0024                                  // "source",
0025                                  "LHC",
0026                                  // "HLT",
0027                                  collectionName,
0028                                  friendlyName,
0029                                  "",
0030                                  edm::TypeWithDict(rawDataType.typeInfo()),
0031                                  false);
0032     desc.setIsProvenanceSetOnRead();
0033     return desc;
0034   }
0035 }  // namespace
0036 
0037 namespace edm {
0038   DaqProvenanceHelper::DaqProvenanceHelper(TypeID const& rawDataType,
0039                                            std::string const& collectionName,
0040                                            std::string const& friendlyName,
0041                                            std::string const& sourceLabel)
0042       : constProductDescription_(
0043             makeDescriptionForDaqProvHelper(rawDataType, collectionName, friendlyName, sourceLabel)),
0044         dummyProvenance_(constProductDescription_.branchID()),
0045         processParameterSet_(),
0046         oldProcessName_(),
0047         oldBranchID_(),
0048         newBranchID_(),
0049         oldProcessHistoryID_(nullptr),
0050         phidMap_() {
0051     // Now we create a process parameter set for the "LHC" process.
0052     // We don't currently use the untracked parameters, However, we make them available, just in case.
0053     std::string const& moduleLabel = constProductDescription_.moduleLabel();
0054     std::string const& processName = constProductDescription_.processName();
0055     typedef std::vector<std::string> vstring;
0056     vstring empty;
0057 
0058     vstring modlbl;
0059     modlbl.reserve(1);
0060     modlbl.push_back(moduleLabel);
0061     processParameterSet_.addParameter("@all_sources", modlbl);
0062 
0063     ParameterSet triggerPaths;
0064     triggerPaths.addParameter<vstring>("@trigger_paths", empty);
0065     processParameterSet_.addParameter<ParameterSet>("@trigger_paths", triggerPaths);
0066 
0067     ParameterSet pseudoInput;
0068     pseudoInput.addParameter<std::string>("@module_edm_type", "Source");
0069     pseudoInput.addParameter<std::string>("@module_label", moduleLabel);
0070     pseudoInput.addParameter<std::string>("@module_type", sourceLabel);
0071     processParameterSet_.addParameter<ParameterSet>(moduleLabel, pseudoInput);
0072 
0073     processParameterSet_.addParameter<vstring>("@all_esmodules", empty);
0074     processParameterSet_.addParameter<vstring>("@all_esprefers", empty);
0075     processParameterSet_.addParameter<vstring>("@all_essources", empty);
0076     processParameterSet_.addParameter<vstring>("@all_loopers", empty);
0077     processParameterSet_.addParameter<vstring>("@all_modules", empty);
0078     processParameterSet_.addParameter<vstring>("@end_paths", empty);
0079     processParameterSet_.addParameter<vstring>("@paths", empty);
0080     processParameterSet_.addParameter<std::string>("@process_name", processName);
0081     // Now we register the process parameter set.
0082     processParameterSet_.registerIt();
0083 
0084     //std::cerr << processParameterSet_.dump() << std::endl;
0085   }
0086 
0087   //default
0088   DaqProvenanceHelper::DaqProvenanceHelper(TypeID const& rawDataType)
0089       : DaqProvenanceHelper(rawDataType, "FEDRawDataCollection", "FEDRawDataCollection", "FedRawDataInputSource") {}
0090 
0091   ProcessHistoryID DaqProvenanceHelper::daqInit(ProductRegistry& productRegistry,
0092                                                 ProcessHistoryRegistry& processHistoryRegistry) const {
0093     // Now we need to set all the metadata
0094     // Add the product to the product registry
0095     productRegistry.copyProduct(constProductDescription_);
0096 
0097     // Insert an entry for this process in the process history registry
0098     // This process is about the data from LHC, and has thus no
0099     // well-defined hardware resources
0100     ProcessHistory ph;
0101     ph.emplace_back(constProductDescription_.processName(),
0102                     processParameterSet_.id(),
0103                     getReleaseVersion(),
0104                     HardwareResourcesDescription());
0105     processHistoryRegistry.registerProcessHistory(ph);
0106 
0107     // Save the process history ID for use every event.
0108     return ph.setProcessHistoryID();
0109   }
0110 
0111   bool DaqProvenanceHelper::matchProcesses(ProcessConfiguration const& newPC, ProcessHistory const& ph) const {
0112     for (auto const& pc : ph) {
0113       if (pc.processName() == oldProcessName_) {
0114         // Assuming here the node hardware information does not matter
0115         return (pc.releaseVersion() == newPC.releaseVersion());
0116       }
0117     }
0118     return false;
0119   }
0120 
0121   void DaqProvenanceHelper::fixMetaData(std::vector<ProcessConfiguration>& pcv, std::vector<ProcessHistory>& phv) {
0122     phv.push_back(ProcessHistory());  // For new processHistory, containing only processConfiguration_
0123     std::vector<ProcessConfiguration> newPCs;
0124     for (auto const& pc : pcv) {
0125       if (pc.processName() == oldProcessName_) {
0126         newPCs.emplace_back(constProductDescription_.processName(),
0127                             processParameterSet_.id(),
0128                             pc.releaseVersion(),
0129                             pc.hardwareResourcesDescription());
0130       }
0131     }
0132     if (newPCs.empty()) {
0133       throw Exception(errors::LogicError) << "\nFatal error in RootFile constructor. Most likely this is because\n"
0134                                           << "the input file contains a FEDRawDataCollection with module label\n"
0135                                           << "\"source\". This is against CMS naming conventions.\n"
0136                                           << "See GitHub Issue 45137 for related details.\n";
0137     }
0138     pcv.reserve(pcv.size() + newPCs.size());
0139     pcv.insert(pcv.end(), newPCs.begin(), newPCs.end());
0140     // update existing process histories
0141     for (auto& ph : phv) {
0142       for (auto const& newPC : newPCs) {
0143         if (ph.empty() || matchProcesses(newPC, ph)) {
0144           ProcessHistoryID oldPHID = ph.id();
0145           ph.push_front(newPC);
0146           ProcessHistoryID newPHID = ph.id();
0147           phidMap_.insert(std::make_pair(oldPHID, newPHID));
0148           break;
0149         }
0150       }
0151     }
0152     // For new process histories, containing only the new process configurations
0153     phv.reserve(phv.size() + newPCs.size());
0154     for (auto const& newPC : newPCs) {
0155       phv.emplace_back();
0156       phv.back().push_front(newPC);
0157     }
0158   }
0159 
0160   void DaqProvenanceHelper::fixMetaData(std::vector<BranchID>& branchID) const {
0161     std::replace(branchID.begin(), branchID.end(), oldBranchID_, newBranchID_);
0162   }
0163 
0164   void DaqProvenanceHelper::fixMetaData(BranchIDLists& branchIDLists) const {
0165     BranchID::value_type oldID = oldBranchID_.id();
0166     BranchID::value_type newID = newBranchID_.id();
0167     for (auto& list : branchIDLists) {
0168       std::replace(list.begin(), list.end(), oldID, newID);
0169     }
0170   }
0171 
0172   void DaqProvenanceHelper::fixMetaData(ProductDependencies& productDependencies) const {
0173     typedef std::map<BranchID, std::set<BranchID> > BCMap;
0174     BCMap& childLookup = productDependencies.mutableChildLookup();
0175     // First fix any old branchID's in the key.
0176     {
0177       BCMap::iterator i = childLookup.find(oldBranchID_);
0178       if (i != childLookup.end()) {
0179         childLookup.insert(std::make_pair(newBranchID_, i->second));
0180         childLookup.erase(i);
0181       }
0182     }
0183     // Now fix any old branchID's in the sets;
0184     for (auto& child : childLookup) {
0185       if (child.second.erase(oldBranchID_) != 0) {
0186         child.second.insert(newBranchID_);
0187       }
0188     }
0189   }
0190 
0191   // Replace process history ID.
0192   ProcessHistoryID const& DaqProvenanceHelper::mapProcessHistoryID(ProcessHistoryID const& phid) {
0193     ProcessHistoryIDMap::const_iterator it = phidMap_.find(phid);
0194     assert(it != phidMap_.end());
0195     oldProcessHistoryID_ = &it->first;
0196     return it->second;
0197   }
0198 
0199   // Replace parentage ID.
0200   ParentageID const& DaqProvenanceHelper::mapParentageID(ParentageID const& parentageID) const {
0201     ParentageIDMap::const_iterator it = parentageIDMap_.find(parentageID);
0202     if (it == parentageIDMap_.end()) {
0203       return parentageID;
0204     }
0205     return it->second;
0206   }
0207 
0208   // Replace branch ID if necessary.
0209   BranchID const& DaqProvenanceHelper::mapBranchID(BranchID const& branchID) const {
0210     return (branchID == oldBranchID_ ? newBranchID_ : branchID);
0211   }
0212 
0213   void DaqProvenanceHelper::setOldParentageIDToNew(ParentageID const& iOld, ParentageID const& iNew) {
0214     parentageIDMap_.insert(std::make_pair(iOld, iNew));
0215   }
0216 
0217 }  // namespace edm