File indexing completed on 2023-03-17 11:03:47
0001 #include <algorithm>
0002 #include <cassert>
0003 #include <vector>
0004
0005 #include "FWCore/Sources/interface/DaqProvenanceHelper.h"
0006
0007 #include "DataFormats/Provenance/interface/BranchChildren.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/GetPassID.h"
0014 #include "FWCore/Reflection/interface/TypeWithDict.h"
0015 #include "FWCore/Version/interface/GetReleaseVersion.h"
0016
0017 namespace {
0018 edm::BranchDescription makeDescriptionForDaqProvHelper(edm::TypeID const& rawDataType,
0019 std::string const& collectionName,
0020 std::string const& friendlyName,
0021 std::string const& sourceLabel) {
0022 edm::BranchDescription desc(edm::InEvent,
0023 "rawDataCollector",
0024
0025 "LHC",
0026
0027 collectionName,
0028 friendlyName,
0029 "",
0030 sourceLabel,
0031 edm::ParameterSetID(),
0032 edm::TypeWithDict(rawDataType.typeInfo()),
0033 false);
0034 desc.setIsProvenanceSetOnRead();
0035 return desc;
0036 }
0037 }
0038
0039 namespace edm {
0040 DaqProvenanceHelper::DaqProvenanceHelper(TypeID const& rawDataType,
0041 std::string const& collectionName,
0042 std::string const& friendlyName,
0043 std::string const& sourceLabel)
0044 : constBranchDescription_(
0045 makeDescriptionForDaqProvHelper(rawDataType, collectionName, friendlyName, sourceLabel)),
0046 dummyProvenance_(constBranchDescription_.branchID()),
0047 processParameterSet_(),
0048 oldProcessName_(),
0049 oldBranchID_(),
0050 newBranchID_(),
0051 oldProcessHistoryID_(nullptr),
0052 phidMap_() {
0053
0054
0055 std::string const& moduleLabel = constBranchDescription_.moduleLabel();
0056 std::string const& processName = constBranchDescription_.processName();
0057 std::string const& moduleName = constBranchDescription_.moduleName();
0058 typedef std::vector<std::string> vstring;
0059 vstring empty;
0060
0061 vstring modlbl;
0062 modlbl.reserve(1);
0063 modlbl.push_back(moduleLabel);
0064 processParameterSet_.addParameter("@all_sources", modlbl);
0065
0066 ParameterSet triggerPaths;
0067 triggerPaths.addParameter<vstring>("@trigger_paths", empty);
0068 processParameterSet_.addParameter<ParameterSet>("@trigger_paths", triggerPaths);
0069
0070 ParameterSet pseudoInput;
0071 pseudoInput.addParameter<std::string>("@module_edm_type", "Source");
0072 pseudoInput.addParameter<std::string>("@module_label", moduleLabel);
0073 pseudoInput.addParameter<std::string>("@module_type", moduleName);
0074 processParameterSet_.addParameter<ParameterSet>(moduleLabel, pseudoInput);
0075
0076 processParameterSet_.addParameter<vstring>("@all_esmodules", empty);
0077 processParameterSet_.addParameter<vstring>("@all_esprefers", empty);
0078 processParameterSet_.addParameter<vstring>("@all_essources", empty);
0079 processParameterSet_.addParameter<vstring>("@all_loopers", empty);
0080 processParameterSet_.addParameter<vstring>("@all_modules", empty);
0081 processParameterSet_.addParameter<vstring>("@end_paths", empty);
0082 processParameterSet_.addParameter<vstring>("@paths", empty);
0083 processParameterSet_.addParameter<std::string>("@process_name", processName);
0084
0085 processParameterSet_.registerIt();
0086
0087
0088 }
0089
0090
0091 DaqProvenanceHelper::DaqProvenanceHelper(TypeID const& rawDataType)
0092 : DaqProvenanceHelper(rawDataType, "FEDRawDataCollection", "FEDRawDataCollection", "FedRawDataInputSource") {}
0093
0094 ProcessHistoryID DaqProvenanceHelper::daqInit(ProductRegistry& productRegistry,
0095 ProcessHistoryRegistry& processHistoryRegistry) const {
0096
0097
0098 productRegistry.copyProduct(constBranchDescription_);
0099
0100
0101 ProcessHistory ph;
0102 ph.emplace_back(constBranchDescription_.processName(), processParameterSet_.id(), getReleaseVersion(), getPassID());
0103 processHistoryRegistry.registerProcessHistory(ph);
0104
0105
0106 return ph.setProcessHistoryID();
0107 }
0108
0109 bool DaqProvenanceHelper::matchProcesses(ProcessConfiguration const& newPC, ProcessHistory const& ph) const {
0110 for (auto const& pc : ph) {
0111 if (pc.processName() == oldProcessName_) {
0112 return (pc.releaseVersion() == newPC.releaseVersion() && pc.passID() == newPC.passID());
0113 }
0114 }
0115 return false;
0116 }
0117
0118 void DaqProvenanceHelper::fixMetaData(std::vector<ProcessConfiguration>& pcv, std::vector<ProcessHistory>& phv) {
0119 phv.push_back(ProcessHistory());
0120 std::vector<ProcessConfiguration> newPCs;
0121 for (auto const& pc : pcv) {
0122 if (pc.processName() == oldProcessName_) {
0123 newPCs.emplace_back(
0124 constBranchDescription_.processName(), processParameterSet_.id(), pc.releaseVersion(), pc.passID());
0125 }
0126 }
0127 assert(!newPCs.empty());
0128 pcv.reserve(pcv.size() + newPCs.size());
0129 pcv.insert(pcv.end(), newPCs.begin(), newPCs.end());
0130
0131 for (auto& ph : phv) {
0132 for (auto const& newPC : newPCs) {
0133 if (ph.empty() || matchProcesses(newPC, ph)) {
0134 ProcessHistoryID oldPHID = ph.id();
0135 ph.push_front(newPC);
0136 ProcessHistoryID newPHID = ph.id();
0137 phidMap_.insert(std::make_pair(oldPHID, newPHID));
0138 break;
0139 }
0140 }
0141 }
0142
0143 phv.reserve(phv.size() + newPCs.size());
0144 for (auto const& newPC : newPCs) {
0145 phv.emplace_back();
0146 phv.back().push_front(newPC);
0147 }
0148 }
0149
0150 void DaqProvenanceHelper::fixMetaData(std::vector<BranchID>& branchID) const {
0151 std::replace(branchID.begin(), branchID.end(), oldBranchID_, newBranchID_);
0152 }
0153
0154 void DaqProvenanceHelper::fixMetaData(BranchIDLists const& branchIDLists) const {
0155 BranchID::value_type oldID = oldBranchID_.id();
0156 BranchID::value_type newID = newBranchID_.id();
0157
0158 BranchIDLists& lists = const_cast<BranchIDLists&>(branchIDLists);
0159 for (auto& list : lists) {
0160 std::replace(list.begin(), list.end(), oldID, newID);
0161 }
0162 }
0163
0164 void DaqProvenanceHelper::fixMetaData(BranchChildren& branchChildren) const {
0165 typedef std::map<BranchID, std::set<BranchID> > BCMap;
0166
0167 BCMap& childLookup = const_cast<BCMap&>(branchChildren.childLookup());
0168
0169 {
0170 BCMap::iterator i = childLookup.find(oldBranchID_);
0171 if (i != childLookup.end()) {
0172 childLookup.insert(std::make_pair(newBranchID_, i->second));
0173 childLookup.erase(i);
0174 }
0175 }
0176
0177 for (auto& child : childLookup) {
0178 if (child.second.erase(oldBranchID_) != 0) {
0179 child.second.insert(newBranchID_);
0180 }
0181 }
0182 }
0183
0184
0185 ProcessHistoryID const& DaqProvenanceHelper::mapProcessHistoryID(ProcessHistoryID const& phid) {
0186 ProcessHistoryIDMap::const_iterator it = phidMap_.find(phid);
0187 assert(it != phidMap_.end());
0188 oldProcessHistoryID_ = &it->first;
0189 return it->second;
0190 }
0191
0192
0193 ParentageID const& DaqProvenanceHelper::mapParentageID(ParentageID const& parentageID) const {
0194 ParentageIDMap::const_iterator it = parentageIDMap_.find(parentageID);
0195 if (it == parentageIDMap_.end()) {
0196 return parentageID;
0197 }
0198 return it->second;
0199 }
0200
0201
0202 BranchID const& DaqProvenanceHelper::mapBranchID(BranchID const& branchID) const {
0203 return (branchID == oldBranchID_ ? newBranchID_ : branchID);
0204 }
0205
0206 void DaqProvenanceHelper::setOldParentageIDToNew(ParentageID const& iOld, ParentageID const& iNew) {
0207 parentageIDMap_.insert(std::make_pair(iOld, iNew));
0208 }
0209
0210 }