CMS 3D CMS Logo

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