CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/FWCore/Sources/src/DaqProvenanceHelper.cc

Go to the documentation of this file.
00001 #include <algorithm>
00002 #include <vector>
00003 
00004 #include "FWCore/Sources/interface/DaqProvenanceHelper.h"
00005 
00006 #include "DataFormats/Provenance/interface/BranchChildren.h"
00007 #include "DataFormats/Provenance/interface/BranchIDList.h"
00008 #include "DataFormats/Provenance/interface/ProcessConfigurationRegistry.h"
00009 #include "DataFormats/Provenance/interface/ProcessHistory.h"
00010 #include "DataFormats/Provenance/interface/ProcessHistoryRegistry.h"
00011 #include "DataFormats/Provenance/interface/ProductRegistry.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/Utilities/interface/GetPassID.h"
00014 #include "FWCore/Version/interface/GetReleaseVersion.h"
00015 
00016 namespace edm {
00017   DaqProvenanceHelper::DaqProvenanceHelper(TypeID const& rawDataType)
00018         : constBranchDescription_(BranchDescription(InEvent
00019                                                   , "rawDataCollector"
00020                                                   //, "source"
00021                                                   , "LHC"
00022                                                   // , "HLT"
00023                                                   , "FEDRawDataCollection"
00024                                                   , "FEDRawDataCollection"
00025                                                   , ""
00026                                                   , "DaqSource"
00027                                                   , ParameterSetID()
00028                                                   , rawDataType
00029                                                   , false))
00030         , dummyProvenance_(constBranchDescription_.branchID())
00031         , processParameterSet_()
00032         , oldProcessName_()
00033         , oldBranchID_()
00034         , newBranchID_()
00035         , oldProcessHistoryID_(0)
00036         , processConfiguration_()
00037         , phidMap_() {
00038     
00039     // Now we create a process parameter set for the "LHC" process.
00040     // We don't currently use the untracked parameters, However, we make them available, just in case.
00041     std::string const& moduleLabel = constBranchDescription_.moduleLabel();
00042     std::string const& processName = constBranchDescription_.processName();
00043     std::string const& moduleName = constBranchDescription_.moduleName();
00044     typedef std::vector<std::string> vstring;
00045     vstring empty;
00046 
00047     vstring modlbl;
00048     modlbl.reserve(1);
00049     modlbl.push_back(moduleLabel);
00050     processParameterSet_.addParameter("@all_sources", modlbl);
00051 
00052     ParameterSet triggerPaths;
00053     triggerPaths.addParameter<vstring>("@trigger_paths", empty);
00054     processParameterSet_.addParameter<ParameterSet>("@trigger_paths", triggerPaths);
00055 
00056     ParameterSet pseudoInput;
00057     pseudoInput.addParameter<std::string>("@module_edm_type", "Source");
00058     pseudoInput.addParameter<std::string>("@module_label", moduleLabel);
00059     pseudoInput.addParameter<std::string>("@module_type", moduleName);
00060     processParameterSet_.addParameter<ParameterSet>(moduleLabel, pseudoInput);
00061 
00062     processParameterSet_.addParameter<vstring>("@all_esmodules", empty);
00063     processParameterSet_.addParameter<vstring>("@all_esprefers", empty);
00064     processParameterSet_.addParameter<vstring>("@all_essources", empty);
00065     processParameterSet_.addParameter<vstring>("@all_loopers", empty);
00066     processParameterSet_.addParameter<vstring>("@all_modules", empty);
00067     processParameterSet_.addParameter<vstring>("@end_paths", empty);
00068     processParameterSet_.addParameter<vstring>("@paths", empty);
00069     processParameterSet_.addParameter<std::string>("@process_name", processName);
00070     // Now we register the process parameter set.
00071     processParameterSet_.registerIt();
00072 
00073     //std::cerr << processParameterSet_.dump() << std::endl;
00074   }
00075 
00076   ProcessHistoryID
00077   DaqProvenanceHelper::daqInit(ProductRegistry& productRegistry) const {
00078     // Now we need to set all the metadata
00079     // Add the product to the product registry  
00080     productRegistry.copyProduct(constBranchDescription_.me());
00081 
00082     // Insert an entry for this process in the process configuration registry
00083     ProcessConfiguration pc(constBranchDescription_.processName(), processParameterSet_.id(), getReleaseVersion(), getPassID());
00084     ProcessConfigurationRegistry::instance()->insertMapped(pc);
00085 
00086     // Insert an entry for this process in the process history registry
00087     ProcessHistory ph;
00088     ph.push_back(pc);
00089     ProcessHistoryRegistry::instance()->insertMapped(ph);
00090 
00091     // Save the process history ID for use every event.
00092     return ph.id();
00093   }
00094 
00095   void
00096   DaqProvenanceHelper::fixMetaData(std::vector<ProcessConfiguration>& pcv) {
00097     bool found = false;
00098     for(std::vector<ProcessConfiguration>::const_iterator it = pcv.begin(), itEnd = pcv.end(); it != itEnd; ++it) {
00099        if(it->processName() == oldProcessName_) {
00100          processConfiguration_ = ProcessConfiguration(constBranchDescription_.processName(),
00101                                                       processParameterSet_.id(),
00102                                                       it->releaseVersion(), it->passID());
00103          pcv.push_back(processConfiguration_);
00104          found = true;
00105          break;
00106        }
00107     }
00108     assert(found);
00109   }
00110 
00111   void
00112   DaqProvenanceHelper::fixMetaData(std::vector<ProcessHistory>& phv) {
00113     phv.reserve(phv.size() + 1);
00114     phv.push_back(ProcessHistory()); // For new processHistory, containing only processConfiguration_.
00115     for(std::vector<ProcessHistory>::iterator it = phv.begin(), itEnd = phv.end(); it != itEnd; ++it) {
00116       ProcessHistoryID oldPHID = it->id();
00117       it->push_front(processConfiguration_);
00118       ProcessHistoryID newPHID = it->id();
00119       phidMap_.insert(std::make_pair(oldPHID, newPHID));
00120     }
00121   }
00122 
00123   void
00124   DaqProvenanceHelper::fixMetaData(std::vector<BranchID>& branchID) const {
00125     std::replace(branchID.begin(), branchID.end(), oldBranchID_, newBranchID_);
00126   }
00127 
00128   void
00129   DaqProvenanceHelper::fixMetaData(BranchIDLists const& branchIDLists) const {
00130     BranchID::value_type oldID = oldBranchID_.id();
00131     BranchID::value_type newID = newBranchID_.id();
00132     // The const_cast is ugly, but it beats the alternatives.
00133     BranchIDLists& lists = const_cast<BranchIDLists&>(branchIDLists);
00134     for(BranchIDLists::iterator it = lists.begin(), itEnd = lists.end(); it != itEnd; ++it) {
00135       std::replace(it->begin(), it->end(), oldID, newID);
00136     }
00137   }
00138 
00139   void
00140   DaqProvenanceHelper::fixMetaData(BranchChildren& branchChildren) const {
00141     typedef std::map<BranchID, std::set<BranchID> > BCMap;
00142     // The const_cast is ugly, but it beats the alternatives.
00143     BCMap& childLookup = const_cast<BCMap&>(branchChildren.childLookup());
00144     // First fix any old branchID's in the key.
00145     {
00146       BCMap::iterator i = childLookup.find(oldBranchID_);
00147       if(i != childLookup.end()) {
00148         childLookup.insert(std::make_pair(newBranchID_, i->second));
00149         childLookup.erase(i);
00150       }
00151     }
00152     // Now fix any old branchID's in the sets;
00153     for(BCMap::iterator it = childLookup.begin(), itEnd = childLookup.end(); it != itEnd; ++it) {
00154       if(it->second.erase(oldBranchID_) != 0) {
00155         it->second.insert(newBranchID_);
00156       }
00157     }
00158   }
00159 
00160   // Replace process history ID.
00161   ProcessHistoryID const&
00162   DaqProvenanceHelper::mapProcessHistoryID(ProcessHistoryID const& phid) {
00163     ProcessHistoryIDMap::const_iterator it = phidMap_.find(phid);
00164     assert(it != phidMap_.end());
00165     oldProcessHistoryID_ = &it->first;
00166     return it->second;
00167   }
00168 
00169   // Replace parentage ID.
00170   ParentageID const&
00171   DaqProvenanceHelper::mapParentageID(ParentageID const& parentageID) const {
00172     ParentageIDMap::const_iterator it = parentageIDMap_.find(parentageID);
00173     if(it == parentageIDMap_.end()) {
00174       return parentageID;
00175     }
00176     return it->second;
00177   }
00178 
00179   // Replace branch ID if necessary.
00180   BranchID const&
00181   DaqProvenanceHelper::mapBranchID(BranchID const& branchID) const {
00182     return(branchID == oldBranchID_ ? newBranchID_ : branchID);
00183   }
00184 }