Go to the documentation of this file.00001 #include "Fireworks/Core/interface/FWLiteJobMetadataManager.h"
00002 #include "Fireworks/Core/interface/FWLiteJobMetadataUpdateRequest.h"
00003 #include "Fireworks/Core/interface/fwLog.h"
00004 #include "Fireworks/Core/interface/FWItemAccessorFactory.h"
00005 #include "DataFormats/Provenance/interface/BranchDescription.h"
00006 #include "DataFormats/FWLite/interface/Event.h"
00007
00008 #include "TFile.h"
00009 #include "TTree.h"
00010 #include <set>
00011
00012 FWLiteJobMetadataManager::FWLiteJobMetadataManager(void)
00013 : FWJobMetadataManager(),
00014 m_event(0)
00015 {}
00016
00022 bool
00023 FWLiteJobMetadataManager::doUpdate(FWJobMetadataUpdateRequest *request)
00024 {
00025 FWLiteJobMetadataUpdateRequest *liteRequest
00026 = dynamic_cast<FWLiteJobMetadataUpdateRequest *>(request);
00027
00028
00029 assert(liteRequest);
00030 if (m_event == liteRequest->event_)
00031 return false;
00032
00033 m_event = liteRequest->event_;
00034 const TFile *file = liteRequest->file_;
00035
00036 assert(file);
00037
00038 usableData().clear();
00039
00040 if (!m_event)
00041 return true;
00042
00043 const std::vector<std::string>& history = m_event->getProcessHistory();
00044
00045
00046
00047 if (history.empty())
00048 fwLog(fwlog::kWarning) << "WARNING: the file '"
00049 << file->GetName() << "' contains no processing history"
00050 " and therefore should have no accessible data";
00051
00052 std::copy(history.rbegin(),history.rend(),
00053 std::back_inserter(processNamesInJob()));
00054
00055 static const std::string s_blank;
00056 const std::vector<edm::BranchDescription>& descriptions =
00057 m_event->getBranchDescriptions();
00058
00059 Data d;
00060
00061
00062 TTree* eventsTree = dynamic_cast<TTree*>(const_cast<TFile*>(file)->Get("Events"));
00063 assert(eventsTree);
00064
00065 std::set<std::string> branchNamesInFile;
00066 TIter nextBranch(eventsTree->GetListOfBranches());
00067 while(TBranch* branch = static_cast<TBranch*>(nextBranch()))
00068 branchNamesInFile.insert(branch->GetName());
00069
00070 typedef std::set<std::string> Purposes;
00071 Purposes purposes;
00072 std::string classType;
00073
00074 for(size_t bi = 0, be = descriptions.size(); bi != be; ++bi)
00075 {
00076 const edm::BranchDescription &desc = descriptions[bi];
00077
00078 if (!desc.present()
00079 || branchNamesInFile.end() == branchNamesInFile.find(desc.branchName()))
00080 continue;
00081
00082 const std::vector<FWRepresentationInfo>& infos
00083 = m_typeAndReps->representationsForType(desc.fullClassName());
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 purposes.clear();
00101 for (size_t ii = 0, ei = infos.size(); ii != ei; ++ii) {
00102
00103
00104
00105 purposes.insert(infos[ii].purpose());
00106 }
00107
00108 if (purposes.empty())
00109 purposes.insert("Table");
00110
00111 for (Purposes::const_iterator itPurpose = purposes.begin(),
00112 itEnd = purposes.end();
00113 itPurpose != itEnd;
00114 ++itPurpose)
00115 {
00116
00117
00118
00119
00120 TClass* theClass = TClass::GetClass(desc.fullClassName().c_str());
00121
00122 if (!theClass)
00123 continue;
00124
00125 if (!theClass->GetTypeInfo())
00126 continue;
00127
00128 const static bool debug = false;
00129
00130 if (!FWItemAccessorFactory::classAccessedAsCollection(theClass) )
00131 {
00132 if (debug) {
00133 fwLog(fwlog::kDebug) << theClass->GetName()
00134 << " will not be displayed in table." << std::endl;
00135 }
00136 continue;
00137 }
00138 d.type_ = desc.fullClassName();
00139 d.purpose_ = *itPurpose;
00140 d.moduleLabel_ = desc.moduleLabel();
00141 d.productInstanceLabel_ = desc.productInstanceName();
00142 d.processName_ = desc.processName();
00143 usableData().push_back(d);
00144 if (debug)
00145 {
00146 fwLog(fwlog::kDebug) << "Add collection will display " << d.type_
00147 << " " << d.moduleLabel_
00148 << " " << d.productInstanceLabel_
00149 << " " << d.processName_ << std::endl;
00150 }
00151 }
00152 }
00153 return true;
00154 }