CMS 3D CMS Logo

RootFile.cc
Go to the documentation of this file.
1 /*----------------------------------------------------------------------
2 ----------------------------------------------------------------------*/
3 
4 #include "RootFile.h"
5 #include "DuplicateChecker.h"
6 #include "InputFile.h"
7 #include "ProvenanceAdaptor.h"
8 #include "RunHelper.h"
9 
45 
48 
49 //used for backward compatibility
54 
55 #include "Rtypes.h"
56 #include "TClass.h"
57 #include "TString.h"
58 #include "TTree.h"
59 #include "TTreeCache.h"
60 
61 #include <algorithm>
62 #include <list>
63 
64 namespace edm {
65 
66  // Algorithm classes for making ProvenanceReader:
68  public:
69  std::unique_ptr<ProvenanceReaderBase> makeReader(RootTree& eventTree, DaqProvenanceHelper const* daqProvenanceHelper) const override;
70  };
72  public:
73  MakeOldProvenanceReader(std::unique_ptr<EntryDescriptionMap>&& entryDescriptionMap) : MakeProvenanceReader(), entryDescriptionMap_(std::move(entryDescriptionMap)) {}
74  std::unique_ptr<ProvenanceReaderBase> makeReader(RootTree& eventTree, DaqProvenanceHelper const* daqProvenanceHelper) const override;
75  private:
77  };
79  public:
80  std::unique_ptr<ProvenanceReaderBase> makeReader(RootTree& eventTree, DaqProvenanceHelper const* daqProvenanceHelper) const override;
81  };
83  public:
84  MakeReducedProvenanceReader(std::vector<ParentageID> const& parentageIDLookup) : parentageIDLookup_(parentageIDLookup) {}
85  std::unique_ptr<ProvenanceReaderBase> makeReader(RootTree& eventTree, DaqProvenanceHelper const* daqProvenanceHelper) const override;
86  private:
87  std::vector<ParentageID> const& parentageIDLookup_;
88  };
89 
90  namespace {
91  void
92  checkReleaseVersion(std::vector<ProcessHistory> processHistoryVector, std::string const& fileName) {
93  std::string releaseVersion = getReleaseVersion();
94  releaseversion::DecomposedReleaseVersion currentRelease(releaseVersion);
95  for(auto const& ph : processHistoryVector) {
96  for(auto const& pc : ph) {
97  if(releaseversion::isEarlierRelease(currentRelease, pc.releaseVersion())) {
99  << "The release you are using, " << getReleaseVersion() << " , predates\n"
100  << "a release (" << pc.releaseVersion() << ") used in writing the input file, " << fileName <<".\n"
101  << "Forward compatibility cannot be supported.\n";
102  }
103  }
104  }
105  }
106  }
107 
108  // This is a helper class for IndexIntoFile.
110  public:
111  explicit RootFileEventFinder(RootTree& eventTree) : eventTree_(eventTree) {}
112  ~RootFileEventFinder() override {}
113 
115  roottree::EntryNumber saveEntry = eventTree_.entryNumber();
116  eventTree_.setEntryNumber(entry);
117  EventAuxiliary eventAux;
118  EventAuxiliary *pEvAux = &eventAux;
119  eventTree_.fillAux<EventAuxiliary>(pEvAux);
120  eventTree_.setEntryNumber(saveEntry);
121  return eventAux.event();
122  }
123 
124  private:
126  };
127 
128 //---------------------------------------------------------------------
130  ProcessConfiguration const& processConfiguration,
131  std::string const& logicalFileName,
132  std::shared_ptr<InputFile> filePtr,
133  std::shared_ptr<EventSkipperByID> eventSkipperByID,
134  bool skipAnyEvents,
135  int remainingEvents,
136  int remainingLumis,
137  unsigned int nStreams,
138  unsigned int treeCacheSize,
139  int treeMaxVirtualSize,
140  InputSource::ProcessingMode processingMode,
141  RunHelperBase* runHelper,
142  bool noEventSort,
143  ProductSelectorRules const& productSelectorRules,
144  InputType inputType,
145  std::shared_ptr<BranchIDListHelper> branchIDListHelper,
146  std::shared_ptr<ThinnedAssociationsHelper> thinnedAssociationsHelper,
147  std::vector<BranchID> const* associationsFromSecondary,
148  std::shared_ptr<DuplicateChecker> duplicateChecker,
149  bool dropDescendants,
150  ProcessHistoryRegistry& processHistoryRegistry,
151  std::vector<std::shared_ptr<IndexIntoFile> > const& indexesIntoFiles,
152  std::vector<std::shared_ptr<IndexIntoFile> >::size_type currentIndexIntoFile,
153  std::vector<ProcessHistoryID>& orderedProcessHistoryIDs,
154  bool bypassVersionCheck,
155  bool labelRawDataLikeMC,
156  bool usingGoToEvent,
157  bool enablePrefetching) :
158  file_(fileName),
159  logicalFile_(logicalFileName),
160  processConfiguration_(processConfiguration),
161  processHistoryRegistry_(&processHistoryRegistry),
162  filePtr_(filePtr),
163  eventSkipperByID_(eventSkipperByID),
164  fileFormatVersion_(),
165  fid_(),
166  indexIntoFileSharedPtr_(new IndexIntoFile),
167  indexIntoFile_(*indexIntoFileSharedPtr_),
168  orderedProcessHistoryIDs_(orderedProcessHistoryIDs),
169  indexIntoFileBegin_(indexIntoFile_.begin(noEventSort ? IndexIntoFile::firstAppearanceOrder : IndexIntoFile::numericalOrder)),
170  indexIntoFileEnd_(indexIntoFileBegin_),
171  indexIntoFileIter_(indexIntoFileBegin_),
172  eventProcessHistoryIDs_(),
173  eventProcessHistoryIter_(eventProcessHistoryIDs_.begin()),
174  savedRunAuxiliary_(),
175  skipAnyEvents_(skipAnyEvents),
176  noEventSort_(noEventSort),
177  whyNotFastClonable_(0),
178  hasNewlyDroppedBranch_(),
179  branchListIndexesUnchanged_(false),
180  eventAux_(),
181  eventTree_(filePtr, InEvent, nStreams, treeMaxVirtualSize, treeCacheSize, roottree::defaultLearningEntries, enablePrefetching, inputType),
182  lumiTree_(filePtr, InLumi, 1, treeMaxVirtualSize, roottree::defaultNonEventCacheSize, roottree::defaultNonEventLearningEntries, enablePrefetching, inputType),
183  runTree_(filePtr, InRun, 1, treeMaxVirtualSize, roottree::defaultNonEventCacheSize, roottree::defaultNonEventLearningEntries, enablePrefetching, inputType),
184  treePointers_(),
185  lastEventEntryNumberRead_(IndexIntoFile::invalidEntry),
186  productRegistry_(),
187  branchIDLists_(),
188  branchIDListHelper_(branchIDListHelper),
189  fileThinnedAssociationsHelper_(),
190  thinnedAssociationsHelper_(thinnedAssociationsHelper),
191  processingMode_(processingMode),
192  runHelper_(runHelper),
193  newBranchToOldBranch_(),
194  eventHistoryTree_(nullptr),
195  eventSelectionIDs_(),
196  branchListIndexes_(),
197  history_(),
198  branchChildren_(new BranchChildren),
199  duplicateChecker_(duplicateChecker),
200  provenanceAdaptor_(),
201  provenanceReaderMaker_(),
202  eventProductProvenanceRetrievers_(),
203  parentageIDLookup_(),
204  daqProvenanceHelper_(),
205  edProductClass_(TClass::GetClass("edm::WrapperBase")) {
206 
207  hasNewlyDroppedBranch_.fill(false);
208 
212 
213  // Read the metadata tree.
214  // We use a smart pointer so the tree will be deleted after use, and not kept for the life of the file.
215  std::unique_ptr<TTree> metaDataTree(dynamic_cast<TTree *>(filePtr_->Get(poolNames::metaDataTreeName().c_str())));
216  if(nullptr == metaDataTree.get()) {
217  throw Exception(errors::FileReadError) << "Could not find tree " << poolNames::metaDataTreeName()
218  << " in the input file.\n";
219  }
220 
221  // To keep things simple, we just read in every possible branch that exists.
222  // We don't pay attention to which branches exist in which file format versions
223 
225  if(metaDataTree->FindBranch(poolNames::fileFormatVersionBranchName().c_str()) != nullptr) {
226  TBranch *fft = metaDataTree->GetBranch(poolNames::fileFormatVersionBranchName().c_str());
227  fft->SetAddress(&fftPtr);
228  roottree::getEntry(fft, 0);
229  metaDataTree->SetBranchAddress(poolNames::fileFormatVersionBranchName().c_str(), &fftPtr);
230  }
231 
232  FileID *fidPtr = &fid_;
233  if(metaDataTree->FindBranch(poolNames::fileIdentifierBranchName().c_str()) != nullptr) {
234  metaDataTree->SetBranchAddress(poolNames::fileIdentifierBranchName().c_str(), &fidPtr);
235  }
236 
237  IndexIntoFile *iifPtr = &indexIntoFile_;
238  if(metaDataTree->FindBranch(poolNames::indexIntoFileBranchName().c_str()) != nullptr) {
239  metaDataTree->SetBranchAddress(poolNames::indexIntoFileBranchName().c_str(), &iifPtr);
240  }
241 
242  // Need to read to a temporary registry so we can do a translation of the BranchKeys.
243  // This preserves backward compatibility against friendly class name algorithm changes.
244  ProductRegistry inputProdDescReg;
245  ProductRegistry *ppReg = &inputProdDescReg;
246  metaDataTree->SetBranchAddress(poolNames::productDescriptionBranchName().c_str(), (&ppReg));
247 
248  typedef std::map<ParameterSetID, ParameterSetBlob> PsetMap;
249  PsetMap psetMap;
250  PsetMap *psetMapPtr = &psetMap;
251  if(metaDataTree->FindBranch(poolNames::parameterSetMapBranchName().c_str()) != nullptr) {
252  //backward compatibility
253  assert(!fileFormatVersion().parameterSetsTree());
254  metaDataTree->SetBranchAddress(poolNames::parameterSetMapBranchName().c_str(), &psetMapPtr);
255  } else {
256  assert(fileFormatVersion().parameterSetsTree());
257  // We use a smart pointer so the tree will be deleted after use, and not kept for the life of the file.
258  std::unique_ptr<TTree> psetTree(dynamic_cast<TTree *>(filePtr_->Get(poolNames::parameterSetsTreeName().c_str())));
259  if(nullptr == psetTree.get()) {
260  throw Exception(errors::FileReadError) << "Could not find tree " << poolNames::parameterSetsTreeName()
261  << " in the input file.\n";
262  }
263 
264  typedef std::pair<ParameterSetID, ParameterSetBlob> IdToBlobs;
265  IdToBlobs idToBlob;
266  IdToBlobs* pIdToBlob = &idToBlob;
267  psetTree->SetBranchAddress(poolNames::idToParameterSetBlobsBranchName().c_str(), &pIdToBlob);
268 
269  std::unique_ptr<TTreeCache> psetTreeCache = roottree::trainCache(psetTree.get(), *filePtr_, roottree::defaultNonEventCacheSize, "*");
270  psetTreeCache->SetEnablePrefetching(false);
271  filePtr_->SetCacheRead(psetTreeCache.get());
272  for(Long64_t i = 0; i != psetTree->GetEntries(); ++i) {
273  psetTree->GetEntry(i);
274  psetMap.insert(idToBlob);
275  }
276  filePtr_->SetCacheRead(nullptr);
277  }
278 
279  // backward compatibility
281  ProcessHistoryRegistry::collection_type *pHistMapPtr = &pHistMap;
282  if(metaDataTree->FindBranch(poolNames::processHistoryMapBranchName().c_str()) != nullptr) {
283  metaDataTree->SetBranchAddress(poolNames::processHistoryMapBranchName().c_str(), &pHistMapPtr);
284  }
285 
287  ProcessHistoryRegistry::vector_type *pHistVectorPtr = &pHistVector;
288  if(metaDataTree->FindBranch(poolNames::processHistoryBranchName().c_str()) != nullptr) {
289  metaDataTree->SetBranchAddress(poolNames::processHistoryBranchName().c_str(), &pHistVectorPtr);
290  }
291 
292  // backward compatibility
293  ProcessConfigurationVector processConfigurations;
294  ProcessConfigurationVector* procConfigVectorPtr = &processConfigurations;
295  if(metaDataTree->FindBranch(poolNames::processConfigurationBranchName().c_str()) != nullptr) {
296  metaDataTree->SetBranchAddress(poolNames::processConfigurationBranchName().c_str(), &procConfigVectorPtr);
297  }
298 
299  auto branchIDListsAPtr = std::make_unique<BranchIDLists>();
300  BranchIDLists* branchIDListsPtr = branchIDListsAPtr.get();
301  if(metaDataTree->FindBranch(poolNames::branchIDListBranchName().c_str()) != nullptr) {
302  metaDataTree->SetBranchAddress(poolNames::branchIDListBranchName().c_str(), &branchIDListsPtr);
303  }
304 
305  if(inputType != InputType::SecondarySource) {
306  fileThinnedAssociationsHelper_ = std::make_unique<ThinnedAssociationsHelper>(); // propagate_const<T> has no reset() function
307  ThinnedAssociationsHelper* thinnedAssociationsHelperPtr = fileThinnedAssociationsHelper_.get();
308  if(metaDataTree->FindBranch(poolNames::thinnedAssociationsHelperBranchName().c_str()) != nullptr) {
309  metaDataTree->SetBranchAddress(poolNames::thinnedAssociationsHelperBranchName().c_str(), &thinnedAssociationsHelperPtr);
310  }
311  }
312 
313  BranchChildren* branchChildrenBuffer = branchChildren_.get();
314  if(metaDataTree->FindBranch(poolNames::productDependenciesBranchName().c_str()) != nullptr) {
315  metaDataTree->SetBranchAddress(poolNames::productDependenciesBranchName().c_str(), &branchChildrenBuffer);
316  }
317 
318  // backward compatibility
319  std::vector<EventProcessHistoryID> *eventHistoryIDsPtr = &eventProcessHistoryIDs_;
320  if(metaDataTree->FindBranch(poolNames::eventHistoryBranchName().c_str()) != nullptr) {
321  metaDataTree->SetBranchAddress(poolNames::eventHistoryBranchName().c_str(), &eventHistoryIDsPtr);
322  }
323 
324  if(metaDataTree->FindBranch(poolNames::moduleDescriptionMapBranchName().c_str()) != nullptr) {
325  if(metaDataTree->GetBranch(poolNames::moduleDescriptionMapBranchName().c_str())->GetSplitLevel() != 0) {
326  metaDataTree->SetBranchStatus((poolNames::moduleDescriptionMapBranchName() + ".*").c_str(), false);
327  } else {
328  metaDataTree->SetBranchStatus(poolNames::moduleDescriptionMapBranchName().c_str(), false);
329  }
330  }
331 
332  // Here we read the metadata tree
333  roottree::getEntry(metaDataTree.get(), 0);
334 
336 
337  // Here we read the event history tree, if we have one.
339 
341  if(!fileFormatVersion().triggerPathsTracked()) {
342  ParameterSetConverter converter(psetMap, psetIdConverter, fileFormatVersion().parameterSetsByReference());
343  } else {
344  // Merge into the parameter set registry.
345  pset::Registry& psetRegistry = *pset::Registry::instance();
346  for(auto const& psetEntry : psetMap) {
347  ParameterSet pset(psetEntry.second.pset());
348  pset.setID(psetEntry.first);
349  // For thread safety, don't update global registries when a secondary source opens a file.
350  if(inputType != InputType::SecondarySource) {
351  psetRegistry.insertMapped(pset);
352  }
353  }
354  }
355  if(!fileFormatVersion().splitProductIDs()) {
356  // Old provenance format input file. Create a provenance adaptor.
357  // propagate_const<T> has no reset() function
358  provenanceAdaptor_ = std::make_unique<ProvenanceAdaptor>(
359  inputProdDescReg, pHistMap, pHistVector, processConfigurations, psetIdConverter, true);
360  // Fill in the branchIDLists branch from the provenance adaptor
361  branchIDLists_ = provenanceAdaptor_->branchIDLists();
362  } else {
363  if(!fileFormatVersion().triggerPathsTracked()) {
364  // New provenance format, but change in ParameterSet Format. Create a provenance adaptor.
365  // propagate_const<T> has no reset() function
366  provenanceAdaptor_ = std::make_unique<ProvenanceAdaptor>(
367  inputProdDescReg, pHistMap, pHistVector, processConfigurations, psetIdConverter, false);
368  }
369  // New provenance format input file. The branchIDLists branch was read directly from the input file.
370  if(metaDataTree->FindBranch(poolNames::branchIDListBranchName().c_str()) == nullptr) {
372  << "Failed to find branchIDLists branch in metaData tree.\n";
373  }
374  branchIDLists_.reset(branchIDListsAPtr.release());
375  }
376 
378  if(metaDataTree->FindBranch(poolNames::thinnedAssociationsHelperBranchName().c_str()) == nullptr) {
380  << "Failed to find thinnedAssociationsHelper branch in metaData tree.\n";
381  }
382  }
383 
384  if(!bypassVersionCheck) {
385  checkReleaseVersion(pHistVector, file());
386  }
387 
388  if(labelRawDataLikeMC) {
389  std::string const rawData("FEDRawDataCollection");
390  std::string const source("source");
391  ProductRegistry::ProductList& pList = inputProdDescReg.productListUpdator();
392  BranchKey finder(rawData, source, "", "");
393  ProductRegistry::ProductList::iterator it = pList.lower_bound(finder);
394  if(it != pList.end() && it->first.friendlyClassName() == rawData && it->first.moduleLabel() == source) {
395  // We found raw data with a module label of source.
396  // We need to change the module label and process name.
397  // Create helper.
398  it->second.init();
399  // propagate_const<T> has no reset() function
400  daqProvenanceHelper_ = std::make_unique<DaqProvenanceHelper>(it->second.unwrappedTypeID());
401  // Create the new branch description
402  BranchDescription const& newBD = daqProvenanceHelper_->branchDescription();
403  // Save info from the old and new branch descriptions
404  daqProvenanceHelper_->saveInfo(it->second, newBD);
405  // Map the new branch name to the old branch name.
406  newBranchToOldBranch_.insert(std::make_pair(newBD.branchName(), it->second.branchName()));
407  // Remove the old branch description from the product Registry.
408  pList.erase(it);
409  // Check that there was only one.
410  it = pList.lower_bound(finder);
411  assert(!(it != pList.end() && it->first.friendlyClassName() == rawData && it->first.moduleLabel() == source));
412  // Insert the new branch description into the product registry.
413  inputProdDescReg.copyProduct(newBD);
414  // Fix up other per file metadata.
415  daqProvenanceHelper_->fixMetaData(processConfigurations, pHistVector);
416  daqProvenanceHelper_->fixMetaData(*branchIDLists_);
417  daqProvenanceHelper_->fixMetaData(*branchChildren_);
418  }
419  }
420 
421  for(auto const& history : pHistVector) {
422  processHistoryRegistry.registerProcessHistory(history);
423  }
424 
426 
427  // Update the branch id info. This has to be done before validateFile since
428  // depending on the file format, the branchIDListHelper_ may have its fixBranchListIndexes call made
429  if(inputType == InputType::Primary) {
431  }
432 
433  validateFile(inputType, usingGoToEvent);
434 
435  // Here, we make the class that will make the ProvenanceReader
436  // It reads whatever trees it needs.
437  // propagate_const<T> has no reset() function
438  provenanceReaderMaker_ = std::unique_ptr<MakeProvenanceReader>(makeProvenanceReaderMaker(inputType).release());
439 
440  // Merge into the hashed registries.
441  if(eventSkipperByID_ && eventSkipperByID_->somethingToSkip()) {
443  }
444 
445  initializeDuplicateChecker(indexesIntoFiles, currentIndexIntoFile);
447  indexIntoFileEnd_ = indexIntoFile_.end(noEventSort ? IndexIntoFile::firstAppearanceOrder : IndexIntoFile::numericalOrder);
449  eventProcessHistoryIter_ = eventProcessHistoryIDs_.begin();
450 
451  // Set product presence information in the product registry.
452  ProductRegistry::ProductList& pList = inputProdDescReg.productListUpdator();
453  for(auto& product : pList) {
454  BranchDescription& prod = product.second;
455  prod.init();
456  treePointers_[prod.branchType()]->setPresence(prod, newBranchToOldBranch(prod.branchName()));
457  }
458 
459  auto newReg = std::make_unique<ProductRegistry>();
460 
461  // Do the translation from the old registry to the new one
462  {
463  ProductRegistry::ProductList const& prodList = inputProdDescReg.productList();
464  for(auto const& product : prodList) {
465  BranchDescription const& prod = product.second;
466  std::string newFriendlyName = friendlyname::friendlyName(prod.className());
467  if(newFriendlyName == prod.friendlyClassName()) {
468  newReg->copyProduct(prod);
469  } else {
470  if(fileFormatVersion().splitProductIDs()) {
472  << "Cannot change friendly class name algorithm without more development work\n"
473  << "to update BranchIDLists and ThinnedAssociationsHelper. Contact the framework group.\n";
474  }
475  BranchDescription newBD(prod);
476  newBD.updateFriendlyClassName();
477  newReg->copyProduct(newBD);
478  newBranchToOldBranch_.insert(std::make_pair(newBD.branchName(), prod.branchName()));
479  }
480  }
481  dropOnInput(*newReg, productSelectorRules, dropDescendants, inputType);
482  if(inputType == InputType::SecondaryFile) {
483  thinnedAssociationsHelper->updateFromSecondaryInput(*fileThinnedAssociationsHelper_,
484  *associationsFromSecondary);
485  } else if (inputType == InputType::Primary) {
486  thinnedAssociationsHelper->updateFromPrimaryInput(*fileThinnedAssociationsHelper_);
487  }
488  // freeze the product registry
489  newReg->setFrozen(inputType != InputType::Primary);
490  productRegistry_.reset(newReg.release());
491  }
492 
493  // Set up information from the product registry.
494  ProductRegistry::ProductList const& prodList = productRegistry()->productList();
495 
496  {
497  std::array<size_t,NumBranchTypes> nBranches;
498  nBranches.fill(0);
499  for(auto const& product : prodList) {
500  ++nBranches[product.second.branchType()];
501  }
502 
503  int i = 0;
504  for(auto t: treePointers_) {
505  t->numberOfBranchesToAdd(nBranches[i]);
506  ++i;
507  }
508  }
509  for(auto const& product : prodList) {
510  BranchDescription const& prod = product.second;
511  treePointers_[prod.branchType()]->addBranch(prod,
513  }
514 
515  // Determine if this file is fast clonable.
516  setIfFastClonable(remainingEvents, remainingLumis);
517 
518  // We are done with our initial reading of EventAuxiliary.
520 
521  // Tell the event tree to begin training at the next read.
523 
524  // Train the run and lumi trees.
525  runTree_.trainCache("*");
526  lumiTree_.trainCache("*");
527  }
528 
530  }
531 
532  void
534  // Called only for old format files.
535  // We use a smart pointer so the tree will be deleted after use, and not kept for the life of the file.
536  std::unique_ptr<TTree> entryDescriptionTree(dynamic_cast<TTree*>(filePtr_->Get(poolNames::entryDescriptionTreeName().c_str())));
537  if(nullptr == entryDescriptionTree.get()) {
538  throw Exception(errors::FileReadError) << "Could not find tree " << poolNames::entryDescriptionTreeName()
539  << " in the input file.\n";
540  }
541 
542  EntryDescriptionID idBuffer;
543  EntryDescriptionID* pidBuffer = &idBuffer;
544  entryDescriptionTree->SetBranchAddress(poolNames::entryDescriptionIDBranchName().c_str(), &pidBuffer);
545 
546  EventEntryDescription entryDescriptionBuffer;
547  EventEntryDescription *pEntryDescriptionBuffer = &entryDescriptionBuffer;
548  entryDescriptionTree->SetBranchAddress(poolNames::entryDescriptionBranchName().c_str(), &pEntryDescriptionBuffer);
549 
550  // Fill in the parentage registry.
552 
553  for(Long64_t i = 0, numEntries = entryDescriptionTree->GetEntries(); i < numEntries; ++i) {
554  roottree::getEntry(entryDescriptionTree.get(), i);
555  if(idBuffer != entryDescriptionBuffer.id()) {
556  throw Exception(errors::EventCorruption) << "Corruption of EntryDescription tree detected.\n";
557  }
558  entryDescriptionMap.insert(std::make_pair(entryDescriptionBuffer.id(),entryDescriptionBuffer));
560  parents.setParents(entryDescriptionBuffer.parents());
562  ParentageID const oldID = parents.id();
563  daqProvenanceHelper_->fixMetaData(parents.parentsForUpdate());
564  ParentageID newID = parents.id();
565  if(newID != oldID) {
566  daqProvenanceHelper_->setOldParentageIDToNew(oldID,newID);
567  }
568  }
569  // For thread safety, don't update global registries when a secondary source opens a file.
570  if(inputType != InputType::SecondarySource) {
571  registry.insertMapped(parents);
572  }
573  }
574  entryDescriptionTree->SetBranchAddress(poolNames::entryDescriptionIDBranchName().c_str(), nullptr);
575  entryDescriptionTree->SetBranchAddress(poolNames::entryDescriptionBranchName().c_str(), nullptr);
576  }
577 
578  void
580  // New format file
581  // We use a smart pointer so the tree will be deleted after use, and not kept for the life of the file.
582  std::unique_ptr<TTree> parentageTree(dynamic_cast<TTree*>(filePtr_->Get(poolNames::parentageTreeName().c_str())));
583  if(nullptr == parentageTree.get()) {
584  throw Exception(errors::FileReadError) << "Could not find tree " << poolNames::parentageTreeName()
585  << " in the input file.\n";
586  }
587 
589  Parentage *pParentageBuffer = &parents;
590  parentageTree->SetBranchAddress(poolNames::parentageBranchName().c_str(), &pParentageBuffer);
591 
593 
594  parentageIDLookup_.reserve(parentageTree->GetEntries());
595  for(Long64_t i = 0, numEntries = parentageTree->GetEntries(); i < numEntries; ++i) {
596  roottree::getEntry(parentageTree.get(), i);
598  ParentageID const oldID = parents.id();
599  daqProvenanceHelper_->fixMetaData(parents.parentsForUpdate());
600  ParentageID newID = parents.id();
601  if(newID != oldID) {
602  daqProvenanceHelper_->setOldParentageIDToNew(oldID,newID);
603  }
604  }
605  // For thread safety, don't update global registries when a secondary source opens a file.
606  if(inputType != InputType::SecondarySource) {
607  registry.insertMapped(parents);
608  }
609  parentageIDLookup_.push_back(parents.id());
610  }
611  parentageTree->SetBranchAddress(poolNames::parentageBranchName().c_str(), nullptr);
612  }
613 
614  void
615  RootFile::setIfFastClonable(int remainingEvents, int remainingLumis) {
616  if(fileFormatVersion().noMetaDataTrees() and !fileFormatVersion().storedProductProvenanceUsed()) {
617  //we must avoid copying the old branch which stored the per product per event provenance
619  return;
620  }
621  if(!fileFormatVersion().splitProductIDs()) {
623  return;
624  }
627  return;
628  }
629  // Find entry for first event in file
631  while(it != indexIntoFileEnd_ && it.getEntryType() != IndexIntoFile::kEvent) {
632  ++it;
633  }
634  if(it == indexIntoFileEnd_) {
636  return;
637  }
638 
639  // From here on, record all reasons we can't fast clone.
643  }
644  if(skipAnyEvents_) {
646  }
647  if(remainingEvents >= 0 && eventTree_.entries() > remainingEvents) {
649  }
650  if(remainingLumis >= 0 && lumiTree_.entries() > remainingLumis) {
652  }
653  if(duplicateChecker_ &&
654  !duplicateChecker_->checkDisabled() &&
655  !duplicateChecker_->noDuplicatesInFile()) {
657  }
658  }
659 
660  std::unique_ptr<FileBlock>
662  return std::make_unique<FileBlock>(fileFormatVersion(),
663  eventTree_.tree(),
665  lumiTree_.tree(),
667  runTree_.tree(),
668  runTree_.metaTree(),
671  file_,
673  modifiedIDs(),
674  branchChildren());
675  }
676 
677  std::string const&
678  RootFile::newBranchToOldBranch(std::string const& newBranch) const {
679  std::map<std::string, std::string>::const_iterator it = newBranchToOldBranch_.find(newBranch);
680  if(it != newBranchToOldBranch_.end()) {
681  return it->second;
682  }
683  return newBranch;
684  }
685 
688  return indexIntoFileIter_;
689  }
690 
691  void
694  }
695 
696  void
697  RootFile::initAssociationsFromSecondary(std::vector<BranchID> const& associationsFromSecondary) {
698  thinnedAssociationsHelper_->initAssociationsFromSecondary(associationsFromSecondary, *fileThinnedAssociationsHelper_);
699  }
700 
701  bool
704  return false;
705  }
706  if(eventSkipperByID_ && eventSkipperByID_->somethingToSkip()) {
707 
708  // See first if the entire lumi or run is skipped, so we won't have to read the event Auxiliary in that case.
710  return true;
711  }
712 
713  // The Lumi is not skipped. If this is an event, see if the event is skipped.
718  eventAux_.id().event())) {
719  return true;
720  }
721  }
722 
723  // Skip runs with no lumis if either lumisToSkip or lumisToProcess have been set to select lumis
725  eventSkipperByID_->skippingLumis()) {
727 
728  // There are no lumis in this run, not even ones we will skip
729  if(iterLumi.peekAheadAtLumi() == IndexIntoFile::invalidLumi) {
730  return true;
731  }
732  // If we get here there are lumis in the run, check to see if we are skipping all of them
733  do {
734  if(!eventSkipperByID_->skipIt(iterLumi.run(), iterLumi.peekAheadAtLumi(), 0U)) {
735  return false;
736  }
737  }
738  while(iterLumi.skipLumiInRun());
739  return true;
740  }
741  }
742  return false;
743  }
744 
745  bool
748  if(duplicateChecker_.get() == nullptr) {
749  return false;
750  }
752  return duplicateChecker_->isDuplicateAndCheckActive(indexIntoFileIter_.processHistoryIDIndex(),
754  }
755 
756  bool
758  return indexIntoFile_.containsItem(run, lumi, event);
759  }
760 
763  // First, account for consecutive skipped entries.
764  while(skipThisEntry()) {
767  }
770  }
771  else {
773  }
774  }
775  // OK, we have an entry that is not skipped.
777  if(entryType == IndexIntoFile::kEnd) {
778  return IndexIntoFile::kEnd;
779  }
780  if(entryType == IndexIntoFile::kRun) {
781  run = indexIntoFileIter_.run();
782  runHelper_->checkForNewRun(run);
783  return IndexIntoFile::kRun;
784  } else if(processingMode_ == InputSource::Runs) {
786  return getNextItemType(run, lumi, event);
787  }
788  if(entryType == IndexIntoFile::kLumi) {
789  run = indexIntoFileIter_.run();
790  lumi = indexIntoFileIter_.lumi();
791  return IndexIntoFile::kLumi;
794  return getNextItemType(run, lumi, event);
795  }
796  if(isDuplicateEvent()) {
798  return getNextItemType(run, lumi, event);
799  }
800  run = indexIntoFileIter_.run();
801  lumi = indexIntoFileIter_.lumi();
803  event = eventAux_.event();
804  return IndexIntoFile::kEvent;
805  }
806 
807  bool
810  itr.advanceToEvent();
811  return itr.getEntryType() == IndexIntoFile::kEnd;
812  }
813 
814  bool
817  int phIndex;
820  IndexIntoFile::EntryNumber_t eventEntry;
821  itr.skipEventBackward(phIndex,
822  run,
823  lumi,
824  eventEntry);
825  itr.skipEventBackward(phIndex,
826  run,
827  lumi,
828  eventEntry);
829  return eventEntry == IndexIntoFile::invalidEntry;
830  }
831 
832  namespace {
833  struct RunItem {
834  RunItem(ProcessHistoryID const& phid, RunNumber_t const& run) :
835  phid_(phid), run_(run) {}
836  ProcessHistoryID phid_;
837  RunNumber_t run_;
838  };
839  struct RunItemSortByRun {
840  bool operator()(RunItem const& a, RunItem const& b) const {
841  return a.run_ < b.run_;
842  }
843  };
844  struct RunItemSortByRunPhid {
845  bool operator()(RunItem const& a, RunItem const& b) const {
846  return a.run_ < b.run_ || (!(b.run_ < a.run_) && a.phid_ < b.phid_);
847  }
848  };
849  struct LumiItem {
850  LumiItem(ProcessHistoryID const& phid, RunNumber_t const& run,
852  phid_(phid), run_(run), lumi_(lumi), firstEventEntry_(entry),
853  lastEventEntry_(entry == IndexIntoFile::invalidEntry ? IndexIntoFile::invalidEntry : entry + 1) {}
854  ProcessHistoryID phid_;
855  RunNumber_t run_;
857  IndexIntoFile::EntryNumber_t firstEventEntry_;
858  IndexIntoFile::EntryNumber_t lastEventEntry_;
859  };
860  struct LumiItemSortByRunLumi {
861  bool operator()(LumiItem const& a, LumiItem const& b) const {
862  return a.run_ < b.run_ || (!(b.run_ < a.run_) && a.lumi_ < b.lumi_);
863  }
864  };
865  struct LumiItemSortByRunLumiPhid {
866  bool operator()(LumiItem const& a, LumiItem const& b) const {
867  if(a.run_ < b.run_) return true;
868  if(b.run_ < a.run_) return false;
869  if(a.lumi_ < b.lumi_) return true;
870  if(b.lumi_ < a.lumi_) return false;
871  return a.phid_ < b.phid_;
872  }
873  };
874  }
875 
876  void
878  // This function is for backward compatibility.
879  // If reading a current format file, indexIntoFile_ is read from the input
880  // file and should always be there. Note that the algorithm below will work
881  // sometimes but often fail with the new format introduced in release 3_8_0.
882  // If it ever becomes necessary to rebuild IndexIntoFile from the new format,
883  // probably a separate function should be written to deal with the task.
884  // This is possible just not implemented yet.
885  assert(!fileFormatVersion().hasIndexIntoFile());
886 
887  typedef std::list<LumiItem> LumiList;
888  LumiList lumis; // (declare 1)
889 
890  typedef std::set<LuminosityBlockID> RunLumiSet;
891  RunLumiSet runLumiSet; // (declare 2)
892 
893  typedef std::list<RunItem> RunList;
894  RunList runs; // (declare 5)
895 
896  typedef std::set<RunNumber_t> RunSet;
897  RunSet runSet; // (declare 4)
898 
899  typedef std::set<RunItem, RunItemSortByRunPhid> RunItemSet;
900  RunItemSet runItemSet; // (declare 3)
901 
902  typedef std::map<RunNumber_t, ProcessHistoryID> PHIDMap;
903  PHIDMap phidMap;
904 
905  RunNumber_t prevRun = 0;
906  LuminosityBlockNumber_t prevLumi = 0;
907  ProcessHistoryID prevPhid;
908  bool iFirst = true;
909 
910  indexIntoFile_.unsortedEventNumbers().clear(); // should already be empty, just being careful
912 
913  // First, loop through the event tree.
914  while(eventTree_.next()) {
915  bool newRun = false;
916  bool newLumi = false;
919 
920  // Save the event numbers as we loop through the event auxiliary to avoid
921  // having to read through the event auxiliary again later. These event numbers
922  // are not actually used in this function, but could be needed elsewhere.
924 
925  ProcessHistoryID reducedPHID = processHistoryRegistry_->reducedProcessHistoryID(eventAux().processHistoryID());
926 
927  if(iFirst || prevPhid != reducedPHID || prevRun != eventAux().run()) {
928  iFirst = false;
929  newRun = newLumi = true;
930  } else if(prevLumi != eventAux().luminosityBlock()) {
931  newLumi = true;
932  }
933  prevPhid = reducedPHID;
934  prevRun = eventAux().run();
935  prevLumi = eventAux().luminosityBlock();
936  if(newLumi) {
937  lumis.emplace_back(reducedPHID,
938  eventAux().run(), eventAux().luminosityBlock(), eventTree_.entryNumber()); // (insert 1)
939  runLumiSet.insert(LuminosityBlockID(eventAux().run(), eventAux().luminosityBlock())); // (insert 2)
940  } else {
941  LumiItem& currentLumi = lumis.back();
942  assert(currentLumi.lastEventEntry_ == eventTree_.entryNumber());
943  ++currentLumi.lastEventEntry_;
944  }
945  if(newRun) {
946  // Insert run in list if it is not already there.
947  RunItem item(reducedPHID, eventAux().run());
948  if(runItemSet.insert(item).second) { // (check 3, insert 3)
949  runs.push_back(std::move(item)); // (insert 5)
950  runSet.insert(eventAux().run()); // (insert 4)
951  phidMap.insert(std::make_pair(eventAux().run(), reducedPHID));
952  }
953  }
954  }
955  // now clean up.
959 
960  // Loop over run entries and fill information.
961 
962  typedef std::map<RunNumber_t, IndexIntoFile::EntryNumber_t> RunMap;
963  RunMap runMap; // (declare 11)
964 
965  typedef std::vector<RunItem> RunVector;
966  RunVector emptyRuns; // (declare 12)
967 
968  if(runTree_.isValid()) {
969  while(runTree_.next()) {
970  // Note: adjacent duplicates will be skipped without an explicit check.
971 
972  std::shared_ptr<RunAuxiliary> runAux = fillRunAuxiliary();
973  ProcessHistoryID reducedPHID = processHistoryRegistry_->reducedProcessHistoryID(runAux->processHistoryID());
974 
975  if(runSet.insert(runAux->run()).second) { // (check 4, insert 4)
976  // This run was not associated with any events.
977  emptyRuns.emplace_back(reducedPHID, runAux->run()); // (insert 12)
978  }
979  runMap.insert(std::make_pair(runAux->run(), runTree_.entryNumber())); // (insert 11)
980  phidMap.insert(std::make_pair(runAux->run(), reducedPHID));
981  }
982  // now clean up.
984  }
985 
986  // Insert the ordered empty runs into the run list.
987  RunItemSortByRun runItemSortByRun;
988  stable_sort_all(emptyRuns, runItemSortByRun);
989 
990  RunList::iterator itRuns = runs.begin(), endRuns = runs.end();
991  for(auto const& emptyRun : emptyRuns) {
992  for(; itRuns != endRuns; ++itRuns) {
993  if(runItemSortByRun(emptyRun, *itRuns)) {
994  break;
995  }
996  }
997  runs.insert(itRuns, emptyRun);
998  }
999 
1000  // Loop over luminosity block entries and fill information.
1001 
1002  typedef std::vector<LumiItem> LumiVector;
1003  LumiVector emptyLumis; // (declare 7)
1004 
1005  typedef std::map<LuminosityBlockID, IndexIntoFile::EntryNumber_t> RunLumiMap;
1006  RunLumiMap runLumiMap; // (declare 6)
1007 
1008  if(lumiTree_.isValid()) {
1009  while(lumiTree_.next()) {
1010  // Note: adjacent duplicates will be skipped without an explicit check.
1011  std::shared_ptr<LuminosityBlockAuxiliary> lumiAux = fillLumiAuxiliary();
1012  LuminosityBlockID lumiID = LuminosityBlockID(lumiAux->run(), lumiAux->luminosityBlock());
1013  if(runLumiSet.insert(lumiID).second) { // (check 2, insert 2)
1014  // This lumi was not associated with any events.
1015  // Use the process history ID from the corresponding run. In cases of practical
1016  // importance, this should be the correct process history ID, but it is possible
1017  // to construct files where this is not the correct process history ID ...
1018  PHIDMap::const_iterator iPhidMap = phidMap.find(lumiAux->run());
1019  assert(iPhidMap != phidMap.end());
1020  emptyLumis.emplace_back(iPhidMap->second, lumiAux->run(), lumiAux->luminosityBlock(), IndexIntoFile::invalidEntry); // (insert 7)
1021  }
1022  runLumiMap.insert(std::make_pair(lumiID, lumiTree_.entryNumber()));
1023  }
1024  // now clean up.
1026  }
1027 
1028  // Insert the ordered empty lumis into the lumi list.
1029  LumiItemSortByRunLumi lumiItemSortByRunLumi;
1030  stable_sort_all(emptyLumis, lumiItemSortByRunLumi);
1031 
1032  LumiList::iterator itLumis = lumis.begin(), endLumis = lumis.end();
1033  for(auto const& emptyLumi : emptyLumis) {
1034  for(; itLumis != endLumis; ++itLumis) {
1035  if(lumiItemSortByRunLumi(emptyLumi, *itLumis)) {
1036  break;
1037  }
1038  }
1039  lumis.insert(itLumis, emptyLumi);
1040  }
1041 
1042  // Create a map of RunItems that gives the order of first appearance in the list.
1043  // Also fill in the vector of process history IDs
1044  typedef std::map<RunItem, int, RunItemSortByRunPhid> RunCountMap;
1045  RunCountMap runCountMap; // Declare (17)
1046  std::vector<ProcessHistoryID>& phids = indexIntoFile_.setProcessHistoryIDs();
1047  assert(phids.empty());
1048  std::vector<IndexIntoFile::RunOrLumiEntry>& entries = indexIntoFile_.setRunOrLumiEntries();
1049  assert(entries.empty());
1050  int rcount = 0;
1051  for(auto& run : runs) {
1052  RunCountMap::const_iterator countMapItem = runCountMap.find(run);
1053  if(countMapItem == runCountMap.end()) {
1054  countMapItem = runCountMap.insert(std::make_pair(run, rcount)).first; // Insert (17)
1055  assert(countMapItem != runCountMap.end());
1056  ++rcount;
1057  }
1058  std::vector<ProcessHistoryID>::const_iterator phidItem = find_in_all(phids, run.phid_);
1059  if(phidItem == phids.end()) {
1060  phids.push_back(run.phid_);
1061  phidItem = phids.end() - 1;
1062  }
1063  entries.emplace_back(
1064  countMapItem->second, // use (17)
1066  runMap[run.run_], // use (11)
1067  phidItem - phids.begin(),
1068  run.run_,
1069  0U,
1072  }
1073 
1074  // Create a map of LumiItems that gives the order of first appearance in the list.
1075  typedef std::map<LumiItem, int, LumiItemSortByRunLumiPhid> LumiCountMap;
1076  LumiCountMap lumiCountMap; // Declare (19)
1077  int lcount = 0;
1078  for(auto& lumi : lumis) {
1079  RunCountMap::const_iterator runCountMapItem = runCountMap.find(RunItem(lumi.phid_, lumi.run_));
1080  assert(runCountMapItem != runCountMap.end());
1081  LumiCountMap::const_iterator countMapItem = lumiCountMap.find(lumi);
1082  if(countMapItem == lumiCountMap.end()) {
1083  countMapItem = lumiCountMap.insert(std::make_pair(lumi, lcount)).first; // Insert (17)
1084  assert(countMapItem != lumiCountMap.end());
1085  ++lcount;
1086  }
1087  std::vector<ProcessHistoryID>::const_iterator phidItem = find_in_all(phids, lumi.phid_);
1088  assert(phidItem != phids.end());
1089  entries.emplace_back(
1090  runCountMapItem->second,
1091  countMapItem->second,
1092  runLumiMap[LuminosityBlockID(lumi.run_, lumi.lumi_)],
1093  phidItem - phids.begin(),
1094  lumi.run_,
1095  lumi.lumi_,
1096  lumi.firstEventEntry_,
1097  lumi.lastEventEntry_);
1098  }
1099  stable_sort_all(entries);
1100  }
1101 
1102  void
1103  RootFile::validateFile(InputType inputType, bool usingGoToEvent) {
1104  if(!fid_.isValid()) {
1106  }
1107  if(!eventTree_.isValid()) {
1109  "'Events' tree is corrupted or not present\n" << "in the input file.\n";
1110  }
1111 
1112  if(fileFormatVersion().hasIndexIntoFile()) {
1113  if(runTree().entries() > 0) {
1114  assert(!indexIntoFile_.empty());
1115  }
1117  if(daqProvenanceHelper_) {
1118  std::vector<ProcessHistoryID>& phidVec = indexIntoFile_.setProcessHistoryIDs();
1119  for(auto& phid : phidVec) {
1120  phid = daqProvenanceHelper_->mapProcessHistoryID(phid);
1121  }
1122  }
1124  }
1125  }
1126  else {
1127  assert(indexIntoFile_.empty());
1129  }
1130 
1133  indexIntoFile_.setEventFinder(std::shared_ptr<IndexIntoFile::EventFinder>(std::make_shared<RootFileEventFinder>(eventTree_)));
1134  // We fill the event numbers explicitly if we need to find events in closed files,
1135  // such as for secondary files (or secondary sources) or if duplicate checking across files.
1136  bool needEventNumbers = false;
1137  bool needIndexesForDuplicateChecker = duplicateChecker_ && duplicateChecker_->checkingAllFiles() && !duplicateChecker_->checkDisabled();
1138  if(inputType != InputType::Primary || needIndexesForDuplicateChecker || usingGoToEvent) {
1139  needEventNumbers = true;
1140  }
1141  bool needEventEntries = false;
1142  if(inputType != InputType::Primary || !noEventSort_) {
1143  // We need event entries for sorting or for secondary files or sources.
1144  needEventEntries = true;
1145  }
1146  indexIntoFile_.fillEventNumbersOrEntries(needEventNumbers, needEventEntries);
1147  }
1148 
1149  void
1151  // Report file opened.
1152  std::string const label = "source";
1153  std::string moduleName = "PoolSource";
1154  filePtr_->inputFileOpened(
1155  logicalFile_,
1156  inputType,
1157  moduleName,
1158  label,
1159  fid_.fid(),
1161  }
1162 
1163  void
1165  // Just to play it safe, zero all pointers to objects in the InputFile to be closed.
1166  eventHistoryTree_ = nullptr;
1167  for(auto& treePointer : treePointers_) {
1168  treePointer->close();
1169  treePointer = nullptr;
1170  }
1171  filePtr_->Close();
1172  filePtr_ = nullptr; // propagate_const<T> has no reset() function
1173  }
1174 
1175  void
1178  // Already read.
1179  return;
1180  }
1182  EventAuxiliary *pEvAux = &eventAux_;
1184  } else {
1185  // for backward compatibility.
1187  EventAux *pEvAux = &eventAux;
1188  eventTree_.fillAux<EventAux>(pEvAux);
1189  conversion(eventAux, eventAux_);
1190  }
1192  }
1193 
1194  bool
1196  if(!eventTree_.current(entry)) {
1197  return false;
1198  }
1199  eventTree_.setEntryNumber(entry);
1201  return true;
1202  }
1203 
1204  void
1206  // We could consider doing delayed reading, but because we have to
1207  // store this History object in a different tree than the event
1208  // data tree, this is too hard to do in this first version.
1209 
1210  if(fileFormatVersion().eventHistoryBranch()) {
1211  // Lumi block number was not in EventID for the relevant releases.
1212  EventID id(eventAux().id().run(), 0, eventAux().id().event());
1213  if(eventProcessHistoryIter_->eventID() != id) {
1216  assert(eventProcessHistoryIter_->eventID() == id);
1217  }
1220  } else if(fileFormatVersion().eventHistoryTree()) {
1221  // for backward compatibility.
1222  History* pHistory = history_.get();
1223  TBranch* eventHistoryBranch = eventHistoryTree_->GetBranch(poolNames::eventHistoryBranchName().c_str());
1224  if(!eventHistoryBranch) {
1226  << "Failed to find history branch in event history tree.\n";
1227  }
1228  eventHistoryBranch->SetAddress(&pHistory);
1230  eventAux_.setProcessHistoryID(history_->processHistoryID());
1231  eventSelectionIDs_.swap(history_->eventSelectionIDs());
1232  branchListIndexes_.swap(history_->branchListIndexes());
1233  } else if(fileFormatVersion().noMetaDataTrees()) {
1234  // Current format
1236  TBranch* eventSelectionIDBranch = eventTree_.tree()->GetBranch(poolNames::eventSelectionsBranchName().c_str());
1237  assert(eventSelectionIDBranch != nullptr);
1238  eventTree_.fillBranchEntry(eventSelectionIDBranch, pESV);
1240  TBranch* branchListIndexesBranch = eventTree_.tree()->GetBranch(poolNames::branchListIndexesBranchName().c_str());
1241  assert(branchListIndexesBranch != nullptr);
1242  eventTree_.fillBranchEntry(branchListIndexesBranch, pBLI);
1243  }
1244  if(provenanceAdaptor_) {
1245  eventAux_.setProcessHistoryID(provenanceAdaptor_->convertID(eventAux().processHistoryID()));
1246  for(auto& esID : eventSelectionIDs_) {
1247  esID = provenanceAdaptor_->convertID(esID);
1248  }
1249  }
1250  if(daqProvenanceHelper_) {
1252  }
1254  // old format. branchListIndexes_ must be filled in from the ProvenanceAdaptor.
1255  provenanceAdaptor_->branchListIndexes(branchListIndexes_);
1256  }
1257  if(branchIDListHelper_) {
1258  branchIDListHelper_->fixBranchListIndexes(branchListIndexes_);
1259  }
1260  }
1261 
1262  std::shared_ptr<LuminosityBlockAuxiliary>
1264  auto lumiAuxiliary = std::make_shared<LuminosityBlockAuxiliary>();
1266  LuminosityBlockAuxiliary *pLumiAux = lumiAuxiliary.get();
1268  } else {
1269  LuminosityBlockAux lumiAux;
1270  LuminosityBlockAux *pLumiAux = &lumiAux;
1272  conversion(lumiAux, *lumiAuxiliary);
1273  }
1274  if(provenanceAdaptor_) {
1275  lumiAuxiliary->setProcessHistoryID(provenanceAdaptor_->convertID(lumiAuxiliary->processHistoryID()));
1276  }
1277  if(daqProvenanceHelper_) {
1278  lumiAuxiliary->setProcessHistoryID(daqProvenanceHelper_->mapProcessHistoryID(lumiAuxiliary->processHistoryID()));
1279  }
1280  if(lumiAuxiliary->luminosityBlock() == 0 && !fileFormatVersion().runsAndLumis()) {
1281  lumiAuxiliary->id() = LuminosityBlockID(RunNumber_t(1), LuminosityBlockNumber_t(1));
1282  }
1283  return lumiAuxiliary;
1284  }
1285 
1286  std::shared_ptr<RunAuxiliary>
1288  auto runAuxiliary = std::make_shared<RunAuxiliary>();
1290  RunAuxiliary *pRunAux = runAuxiliary.get();
1291  runTree_.fillAux<RunAuxiliary>(pRunAux);
1292  } else {
1293  RunAux runAux;
1294  RunAux *pRunAux = &runAux;
1295  runTree_.fillAux<RunAux>(pRunAux);
1296  conversion(runAux, *runAuxiliary);
1297  }
1298  if(provenanceAdaptor_) {
1299  runAuxiliary->setProcessHistoryID(provenanceAdaptor_->convertID(runAuxiliary->processHistoryID()));
1300  }
1301  if(daqProvenanceHelper_) {
1302  runAuxiliary->setProcessHistoryID(daqProvenanceHelper_->mapProcessHistoryID(runAuxiliary->processHistoryID()));
1303  }
1304  return runAuxiliary;
1305  }
1306 
1307  bool
1309  while(offset > 0 && indexIntoFileIter_ != indexIntoFileEnd_) {
1310 
1311  int phIndexOfSkippedEvent = IndexIntoFile::invalidIndex;
1312  RunNumber_t runOfSkippedEvent = IndexIntoFile::invalidRun;
1315 
1316  indexIntoFileIter_.skipEventForward(phIndexOfSkippedEvent,
1317  runOfSkippedEvent,
1318  lumiOfSkippedEvent,
1319  skippedEventEntry);
1320 
1321  // At the end of the file and there were no more events to skip
1322  if(skippedEventEntry == IndexIntoFile::invalidEntry) break;
1323 
1324  if(eventSkipperByID_ && eventSkipperByID_->somethingToSkip()) {
1325  fillEventAuxiliary(skippedEventEntry);
1326  if(eventSkipperByID_->skipIt(runOfSkippedEvent, lumiOfSkippedEvent, eventAux_.id().event())) {
1327  continue;
1328  }
1329  }
1330  if(duplicateChecker_ &&
1331  !duplicateChecker_->checkDisabled() &&
1332  !duplicateChecker_->noDuplicatesInFile()) {
1333 
1334  fillEventAuxiliary(skippedEventEntry);
1335  if(duplicateChecker_->isDuplicateAndCheckActive(phIndexOfSkippedEvent,
1336  runOfSkippedEvent,
1337  lumiOfSkippedEvent,
1338  eventAux_.id().event(),
1339  file_)) {
1340  continue;
1341  }
1342  }
1343  --offset;
1344  }
1345 
1346  while(offset < 0) {
1347 
1348  if(duplicateChecker_) {
1349  duplicateChecker_->disable();
1350  }
1351 
1352  int phIndexOfEvent = IndexIntoFile::invalidIndex;
1356 
1357  indexIntoFileIter_.skipEventBackward(phIndexOfEvent,
1358  runOfEvent,
1359  lumiOfEvent,
1360  eventEntry);
1361 
1362  if(eventEntry == IndexIntoFile::invalidEntry) break;
1363 
1364  if(eventSkipperByID_ && eventSkipperByID_->somethingToSkip()) {
1365  fillEventAuxiliary(eventEntry);
1366  if(eventSkipperByID_->skipIt(runOfEvent, lumiOfEvent, eventAux_.id().event())) {
1367  continue;
1368  }
1369  }
1370  ++offset;
1371  }
1373  }
1374 
1375  bool
1377 
1379 
1380  if(duplicateChecker_) {
1381  duplicateChecker_->disable();
1382  }
1383 
1386 
1388  indexIntoFile_.findPosition(sortOrder, eventID.run(), eventID.luminosityBlock(), eventID.event());
1389 
1390  if(iter == indexIntoFile_.end(sortOrder)) {
1391  return false;
1392  }
1393  indexIntoFileIter_ = iter;
1394  return true;
1395  }
1396 
1397  // readEvent() is responsible for creating, and setting up, the
1398  // EventPrincipal.
1399  //
1400  // 1. create an EventPrincipal with a unique EventID
1401  // 2. For each entry in the provenance, put in one ProductResolver,
1402  // holding the Provenance for the corresponding EDProduct.
1403  // 3. set up the the EventPrincipal to know about this ProductResolver.
1404  //
1405  // We do *not* create the EDProduct instance (the equivalent of reading
1406  // the branch containing this EDProduct. That will be done by the Delayed Reader,
1407  // when it is asked to do so.
1408  //
1409  void
1413  // read the event auxiliary if not alrady read.
1415 
1416  // read the event
1417  readCurrentEvent(principal);
1418 
1419  runHelper_->checkRunConsistency(eventAux().run(), indexIntoFileIter_.run());
1420  runHelper_->checkLumiConsistency(eventAux().luminosityBlock(), indexIntoFileIter_.lumi());
1421 
1422  // If this next assert shows up in performance profiling or significantly affects memory, then these three lines should be deleted.
1423  // The IndexIntoFile should guarantee that it never fails.
1425  ProcessHistoryID const& reducedPHID = processHistoryRegistry_->reducedProcessHistoryID(idToCheck);
1427 
1429  }
1430 
1431  // Reads event at the current entry in the event tree
1432  bool
1434  if(!eventTree_.current()) {
1435  return false;
1436  }
1438  if(!fileFormatVersion().lumiInEventID()) {
1439  //ugly, but will disappear when the backward compatibility is done with schema evolution.
1440  const_cast<EventID&>(eventAux_.id()).setLuminosityBlockNumber(eventAux_.oldLuminosityBlock());
1442  }
1443  fillEventHistory();
1444  runHelper_->overrideRunNumber(eventAux_.id(), eventAux().isRealData());
1445 
1446  // We're not done ... so prepare the EventPrincipal
1448  principal.fillEventPrincipal(eventAux(),
1452  *(makeProductProvenanceRetriever(principal.streamID().value())),
1454 
1455  // report event read from file
1456  filePtr_->eventReadFromFile();
1457  return true;
1458  }
1459 
1460  void
1462  eventTree_.setEntryNumber(entry);
1463  }
1464 
1465  std::shared_ptr<RunAuxiliary>
1467  if(runHelper_->fakeNewRun()) {
1468  runHelper_->overrideRunNumber(savedRunAuxiliary_->id());
1469  return savedRunAuxiliary();
1470  }
1473 
1474  // Begin code for backward compatibility before the existence of run trees.
1475  if(!runTree_.isValid()) {
1476 
1477  // prior to the support of run trees.
1478  // RunAuxiliary did not contain a valid timestamp. Take it from the next event.
1480  assert(eventEntry != IndexIntoFile::invalidEntry);
1481  assert(eventTree_.current(eventEntry));
1482  fillEventAuxiliary(eventEntry);
1483 
1485  runHelper_->overrideRunNumber(run);
1486  savedRunAuxiliary_ = std::make_shared<RunAuxiliary>(run.run(), eventAux().time(), Timestamp::invalidTimestamp());
1487  return savedRunAuxiliary();
1488  }
1489  // End code for backward compatibility before the existence of run trees.
1491  std::shared_ptr<RunAuxiliary> runAuxiliary = fillRunAuxiliary();
1492  assert(runAuxiliary->run() == indexIntoFileIter_.run());
1493  runHelper_->overrideRunNumber(runAuxiliary->id());
1494  filePtr_->reportInputRunNumber(runAuxiliary->run());
1495  // If RunAuxiliary did not contain a valid begin timestamp, invalidate any end timestamp.
1496  if(runAuxiliary->beginTime() == Timestamp::invalidTimestamp()) {
1497  runAuxiliary->setEndTime(Timestamp::invalidTimestamp());
1498  }
1499 
1500  // If RunAuxiliary did not contain a valid timestamp, or if this an old format file from
1501  // when the Run's ProcessHistory included only processes where products were added to the Run itself,
1502  // we attempt to read the first event in the run to get appropriate info.
1503  if(runAuxiliary->beginTime() == Timestamp::invalidTimestamp() ||
1505 
1507  // If we have a valid event, use its information.
1508  if(eventEntry != IndexIntoFile::invalidEntry) {
1509  assert(eventTree_.current(eventEntry));
1510  fillEventAuxiliary(eventEntry);
1511 
1512  // RunAuxiliary did not contain a valid timestamp. Take it from the next event in this run if there is one.
1513  if(runAuxiliary->beginTime() == Timestamp::invalidTimestamp()) {
1514  runAuxiliary->setBeginTime(eventAux().time());
1515  }
1516 
1517  // For backwards compatibility when the Run's ProcessHistory included only processes where products were added to the
1518  // Run, and then the Run and Event auxiliaries could be different. Use the event ProcessHistoryID if there is one. It should
1519  // almost always be correct by the current definition (processes included if any products are added. This makes the run, lumi,
1520  // and event ProcessHistory's always be the same if no file merging occurs).
1521  if(!fileFormatVersion().processHistorySameWithinRun()) {
1522  fillEventHistory();
1523  runAuxiliary->setProcessHistoryID(eventAux().processHistoryID());
1524  }
1525  }
1526  }
1527  savedRunAuxiliary_ = runAuxiliary;
1528  return runAuxiliary;
1529  }
1530 
1531  void
1533  if(!runHelper_->fakeNewRun()) {
1537  }
1538  // Begin code for backward compatibility before the existence of run trees.
1539  if(!runTree_.isValid()) {
1540  return;
1541  }
1542  // End code for backward compatibility before the existence of run trees.
1543  // NOTE: we use 0 for the index since do not do delayed reads for RunPrincipals
1546  // Read in all the products now.
1547  runPrincipal.readAllFromSourceAndMergeImmediately();
1548  }
1549 
1550 
1551  std::shared_ptr<LuminosityBlockAuxiliary>
1555  // Begin code for backward compatibility before the existence of lumi trees.
1556  if(!lumiTree_.isValid()) {
1558  assert(eventEntry != IndexIntoFile::invalidEntry);
1559  assert(eventTree_.current(eventEntry));
1560  fillEventAuxiliary(eventEntry);
1561 
1563  runHelper_->overrideRunNumber(lumi);
1564  return std::make_shared<LuminosityBlockAuxiliary>(lumi.run(), lumi.luminosityBlock(), eventAux().time(), Timestamp::invalidTimestamp());
1565  }
1566  // End code for backward compatibility before the existence of lumi trees.
1568  std::shared_ptr<LuminosityBlockAuxiliary> lumiAuxiliary = fillLumiAuxiliary();
1569  assert(lumiAuxiliary->run() == indexIntoFileIter_.run());
1570  assert(lumiAuxiliary->luminosityBlock() == indexIntoFileIter_.lumi());
1571  runHelper_->overrideRunNumber(lumiAuxiliary->id());
1572  filePtr_->reportInputLumiSection(lumiAuxiliary->run(), lumiAuxiliary->luminosityBlock());
1573  if(lumiAuxiliary->beginTime() == Timestamp::invalidTimestamp()) {
1575  if(eventEntry != IndexIntoFile::invalidEntry) {
1576  assert(eventTree_.current(eventEntry));
1577  fillEventAuxiliary(eventEntry);
1578 
1579  lumiAuxiliary->setBeginTime(eventAux().time());
1580  }
1581  lumiAuxiliary->setEndTime(Timestamp::invalidTimestamp());
1582  }
1583  if(!fileFormatVersion().processHistorySameWithinRun() && savedRunAuxiliary_) {
1584  lumiAuxiliary->setProcessHistoryID(savedRunAuxiliary_->processHistoryID());
1585  }
1586  return lumiAuxiliary;
1587  }
1588 
1589  void
1593  // Begin code for backward compatibility before the existence of lumi trees.
1594  if(!lumiTree_.isValid()) {
1596  return;
1597  }
1598  // End code for backward compatibility before the existence of lumi trees.
1600  // NOTE: we use 0 for the index since do not do delayed reads for LuminosityBlockPrincipals
1603  // Read in all the products now.
1604  lumiPrincipal.readAllFromSourceAndMergeImmediately();
1606  }
1607 
1608  bool
1611  if(indexIntoFileIter_ == indexIntoFileEnd_) return false;
1613  return true;
1614  }
1615 
1616  bool
1619  if(indexIntoFileIter_ == indexIntoFileEnd_) return false;
1621  return true;
1622  }
1623 
1624  bool
1627  if(indexIntoFileIter_ == indexIntoFileEnd_) return false;
1629  return true;
1630  }
1631 
1632  bool
1636  }
1639  if(run != indexIntoFileIter_.run()) return false;
1640  if(lumi != indexIntoFileIter_.lumi()) return false;
1642  return true;
1643  }
1644 
1645  void
1647  // Read in the event history tree, if we have one...
1648  if(fileFormatVersion().eventHistoryTree()) {
1649  history_ = std::make_unique<History>(); // propagate_const<T> has no reset() function
1650  eventHistoryTree_ = dynamic_cast<TTree*>(filePtr_->Get(poolNames::eventHistoryTreeName().c_str()));
1651  if(!eventHistoryTree_) {
1653  << "Failed to find the event history tree.\n";
1654  }
1655  }
1656  }
1657 
1658  void
1660  std::vector<std::shared_ptr<IndexIntoFile> > const& indexesIntoFiles,
1661  std::vector<std::shared_ptr<IndexIntoFile> >::size_type currentIndexIntoFile) {
1662  if(duplicateChecker_) {
1663  if(eventTree_.next()) {
1665  duplicateChecker_->inputFileOpened(eventAux().isRealData(),
1667  indexesIntoFiles,
1668  currentIndexIntoFile);
1669  }
1671  }
1672  }
1673 
1674  void
1675  RootFile::markBranchToBeDropped(bool dropDescendants, BranchDescription const& branch, std::set<BranchID>& branchesToDrop,
1676  std::map<BranchID, BranchID> const& droppedToKeptAlias) const {
1677  if(dropDescendants) {
1678  branchChildren_->appendToDescendants(branch, branchesToDrop, droppedToKeptAlias);
1679  } else {
1680  branchesToDrop.insert(branch.branchID());
1681  }
1682  }
1683 
1684  void
1685  RootFile::dropOnInput (ProductRegistry& reg, ProductSelectorRules const& rules, bool dropDescendants, InputType inputType) {
1686 
1687  // This is the selector for drop on input.
1688  ProductSelector productSelector;
1689  productSelector.initialize(rules, reg.allBranchDescriptions());
1690 
1691  std::vector<BranchDescription const*> associationDescriptions;
1692 
1694  // Do drop on input. On the first pass, just fill in a set of branches to be dropped.
1695  std::set<BranchID> branchesToDrop;
1696  std::map<BranchID, BranchID> droppedToKeptAlias;
1697  for(auto const& product : prodList) {
1698  BranchDescription const& prod = product.second;
1699  if (prod.branchID() != prod.originalBranchID() && prod.present()) {
1700  droppedToKeptAlias[prod.originalBranchID()] = prod.branchID();
1701  }
1702  }
1703  for(auto const& product : prodList) {
1704  BranchDescription const& prod = product.second;
1705  // Special handling for ThinnedAssociations
1706  if(prod.unwrappedType() == typeid(ThinnedAssociation) && prod.present()) {
1707  if(inputType != InputType::SecondarySource) {
1708  associationDescriptions.push_back(&prod);
1709  } else {
1710  markBranchToBeDropped(dropDescendants, prod, branchesToDrop, droppedToKeptAlias);
1711  }
1712  } else if(!productSelector.selected(prod)) {
1713  markBranchToBeDropped(dropDescendants, prod, branchesToDrop, droppedToKeptAlias);
1714  }
1715  }
1716 
1717  if(inputType != InputType::SecondarySource) {
1718 
1719  // Decide whether to keep the thinned associations and corresponding
1720  // entries in the helper. For secondary source they are all dropped,
1721  // but in other cases we look for thinned collections the associations
1722  // redirect a Ref or Ptr to when dereferencing them.
1723 
1724  // Need a list of kept products in order to determine which thinned associations
1725  // are kept.
1726  std::set<BranchID> keptProductsInEvent;
1727  for(auto const& product : prodList) {
1728  BranchDescription const& prod = product.second;
1729  if( branchesToDrop.find(prod.branchID()) == branchesToDrop.end() &&
1730  prod.present() &&
1731  prod.branchType() == InEvent) {
1732  keptProductsInEvent.insert(prod.branchID());
1733  }
1734  }
1735 
1736  // Decide which ThinnedAssociations to keep and store the decision in keepAssociation
1737  std::map<BranchID, bool> keepAssociation;
1738  fileThinnedAssociationsHelper_->selectAssociationProducts(associationDescriptions,
1739  keptProductsInEvent,
1740  keepAssociation);
1741 
1742  for(auto association : associationDescriptions) {
1743  if(!keepAssociation[association->branchID()]) {
1744  markBranchToBeDropped(dropDescendants, *association, branchesToDrop, droppedToKeptAlias);
1745  }
1746  }
1747 
1748  // Also delete the dropped associations from the ThinnedAssociationsHelper
1749  auto temp = std::make_unique<ThinnedAssociationsHelper>();
1750  for(auto const& associationBranches : fileThinnedAssociationsHelper_->data()) {
1751  auto iter = keepAssociation.find(associationBranches.association());
1752  if(iter != keepAssociation.end() && iter->second) {
1753  temp->addAssociation(associationBranches);
1754  }
1755  }
1756  // propagate_const<T> has no reset() function
1757  fileThinnedAssociationsHelper_ = std::unique_ptr<ThinnedAssociationsHelper>(temp.release());
1758  }
1759 
1760  // On this pass, actually drop the branches.
1761  std::set<BranchID>::const_iterator branchesToDropEnd = branchesToDrop.end();
1762  for(ProductRegistry::ProductList::iterator it = prodList.begin(), itEnd = prodList.end(); it != itEnd;) {
1763  BranchDescription const& prod = it->second;
1764  bool drop = branchesToDrop.find(prod.branchID()) != branchesToDropEnd;
1765  if(drop) {
1766  if(!prod.dropped()) {
1767  if(productSelector.selected(prod) && prod.unwrappedType() != typeid(ThinnedAssociation)) {
1768  LogWarning("RootFile")
1769  << "Branch '" << prod.branchName() << "' is being dropped from the input\n"
1770  << "of file '" << file_ << "' because it is dependent on a branch\n"
1771  << "that was explicitly dropped.\n";
1772  }
1773  treePointers_[prod.branchType()]->dropBranch(newBranchToOldBranch(prod.branchName()));
1774  hasNewlyDroppedBranch_[prod.branchType()] = true;
1775  }
1776  ProductRegistry::ProductList::iterator icopy = it;
1777  ++it;
1778  prodList.erase(icopy);
1779  } else {
1780  ++it;
1781  }
1782  }
1783 
1784  // Drop on input mergeable run and lumi products, this needs to be invoked for secondary file input
1785  if(inputType == InputType::SecondaryFile) {
1786  TString tString;
1787  for(ProductRegistry::ProductList::iterator it = prodList.begin(), itEnd = prodList.end(); it != itEnd;) {
1788  BranchDescription const& prod = it->second;
1789  if(prod.branchType() != InEvent) {
1790  TClass* cp = TClass::GetClass(prod.wrappedName().c_str());
1791  void* p = cp->New();
1792  int offset = cp->GetBaseClassOffset(edProductClass_);
1793  std::unique_ptr<WrapperBase> edp = getWrapperBasePtr(p, offset);
1794  if(edp->isMergeable()) {
1795  treePointers_[prod.branchType()]->dropBranch(newBranchToOldBranch(prod.branchName()));
1796  ProductRegistry::ProductList::iterator icopy = it;
1797  ++it;
1798  prodList.erase(icopy);
1799  } else {
1800  ++it;
1801  }
1802  }
1803  else ++it;
1804  }
1805  }
1806  }
1807 
1808  void
1809  RootFile::setSignals(signalslot::Signal<void(StreamContext const&, ModuleCallingContext const&)> const* preEventReadSource,
1810  signalslot::Signal<void(StreamContext const&, ModuleCallingContext const&)> const* postEventReadSource) {
1811  eventTree_.setSignals(preEventReadSource,postEventReadSource);
1812  }
1813 
1814 
1815  std::unique_ptr<MakeProvenanceReader>
1818  readParentageTree(inputType);
1819  return std::make_unique<MakeReducedProvenanceReader>(parentageIDLookup_);
1820  } else if(fileFormatVersion_.splitProductIDs()) {
1821  readParentageTree(inputType);
1822  return std::make_unique<MakeFullProvenanceReader>();
1823  } else if(fileFormatVersion_.perEventProductIDs()) {
1824  auto entryDescriptionMap = std::make_unique<EntryDescriptionMap>();
1825  readEntryDescriptionTree(*entryDescriptionMap, inputType);
1826  return std::make_unique<MakeOldProvenanceReader>(std::move(entryDescriptionMap));
1827  } else {
1828  return std::make_unique<MakeDummyProvenanceReader>();
1829  }
1830  }
1831 
1832  std::shared_ptr<ProductProvenanceRetriever>
1834  if(eventProductProvenanceRetrievers_.size()<=iStreamID) {
1835  eventProductProvenanceRetrievers_.resize(iStreamID+1);
1836  }
1837  if(!eventProductProvenanceRetrievers_[iStreamID]) {
1838  // propagate_const<T> has no reset() function
1839  eventProductProvenanceRetrievers_[iStreamID] = std::make_shared<ProductProvenanceRetriever>(provenanceReaderMaker_->makeReader(eventTree_, daqProvenanceHelper_.get()));
1840  }
1841  eventProductProvenanceRetrievers_[iStreamID]->reset();
1842  return eventProductProvenanceRetriever(iStreamID);
1843  }
1844 
1846  public:
1847  ReducedProvenanceReader(RootTree* iRootTree, std::vector<ParentageID> const& iParentageIDLookup, DaqProvenanceHelper const* daqProvenanceHelper);
1848 
1849  std::set<ProductProvenance> readProvenance(unsigned int) const override;
1850 private:
1851  void readProvenanceAsync(WaitingTask* task,
1852  ModuleCallingContext const* moduleCallingContext,
1853  unsigned int transitionIndex,
1854  std::atomic<const std::set<ProductProvenance>*>& writeTo
1855  ) const override;
1856 
1861  std::vector<ParentageID> const& parentageIDLookup_;
1863  std::shared_ptr<std::recursive_mutex> mutex_;
1865  };
1866 
1868  RootTree* iRootTree,
1869  std::vector<ParentageID> const& iParentageIDLookup,
1870  DaqProvenanceHelper const* daqProvenanceHelper) :
1872  rootTree_(iRootTree),
1873  pProvVector_(&provVector_),
1874  parentageIDLookup_(iParentageIDLookup),
1875  daqProvenanceHelper_(daqProvenanceHelper),
1876  mutex_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().second),
1877  acquirer_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().first)
1878  {
1879  provBranch_ = rootTree_->tree()->GetBranch(BranchTypeToProductProvenanceBranchName(rootTree_->branchType()).c_str());
1880  }
1881 
1882  namespace {
1884  template<typename R>
1885  void readProvenanceAsyncImpl(R const* iThis,
1887  WaitingTask* task,
1888  unsigned int transitionIndex,
1889  std::atomic<const std::set<ProductProvenance>*>& writeTo,
1890  ModuleCallingContext const* iContext,
1891  SignalType const* pre,
1892  SignalType const* post) {
1893  if(nullptr == writeTo.load()) {
1894  //need to be sure the task isn't run until after the read
1895  WaitingTaskHolder taskHolder{task};
1896  auto pWriteTo = &writeTo;
1897 
1898  auto serviceToken = ServiceRegistry::instance().presentToken();
1899 
1900  chain.push([holder = std::move(taskHolder), pWriteTo,iThis,transitionIndex, iContext, pre,post, serviceToken]() mutable {
1901 
1902  if(nullptr == pWriteTo->load()) {
1903  ServiceRegistry::Operate operate(serviceToken);
1904  std::unique_ptr<const std::set<ProductProvenance>> prov;
1905  try {
1906  if(pre) {pre->emit(*(iContext->getStreamContext()),*iContext);}
1907  prov = std::make_unique<const std::set<ProductProvenance>>(iThis->readProvenance(transitionIndex));
1908  if(post) {post->emit(*(iContext->getStreamContext()),*iContext);}
1909 
1910  }catch(...) {
1911  if(post) {post->emit(*(iContext->getStreamContext()),*iContext);}
1912 
1913  holder.doneWaiting(std::current_exception());
1914  return;
1915  }
1916  const std::set<ProductProvenance>* expected = nullptr;
1917 
1918  if(pWriteTo->compare_exchange_strong(expected,prov.get())) {
1919  prov.release();
1920  }
1921  }
1922  holder.doneWaiting(std::exception_ptr());
1923  });
1924  }
1925  }
1926  }
1927 
1928  void
1930  ModuleCallingContext const* moduleCallingContext,
1931  unsigned int transitionIndex,
1932  std::atomic<const std::set<ProductProvenance>*>& writeTo
1933  ) const {
1934  readProvenanceAsyncImpl(this,
1936  task,
1937  transitionIndex,
1938  writeTo,
1939  moduleCallingContext,
1940  rootTree_->rootDelayedReader()->preEventReadFromSourceSignal(),
1941  rootTree_->rootDelayedReader()->postEventReadFromSourceSignal());
1942  }
1943 
1944  std::set<ProductProvenance>
1945  ReducedProvenanceReader::readProvenance(unsigned int transitionIndex) const {
1946  {
1947  std::lock_guard<std::recursive_mutex> guard(*mutex_);
1948  ReducedProvenanceReader* me = const_cast<ReducedProvenanceReader*>(this);
1949  me->rootTree_->fillBranchEntry(me->provBranch_, me->rootTree_->entryNumberForIndex(transitionIndex), me->pProvVector_);
1950  }
1951  std::set<ProductProvenance> retValue;
1952  if(daqProvenanceHelper_) {
1953  for(auto const& prov : provVector_) {
1954  BranchID bid(prov.branchID_);
1955  retValue.emplace(daqProvenanceHelper_->mapBranchID(BranchID(prov.branchID_)),
1956  daqProvenanceHelper_->mapParentageID(parentageIDLookup_[prov.parentageIDIndex_]));
1957  }
1958  } else {
1959  for(auto const& prov : provVector_) {
1960  if(prov.parentageIDIndex_ >= parentageIDLookup_.size()) {
1962  << "ReducedProvenanceReader::ReadProvenance\n"
1963  << "The parentage ID index value " << prov.parentageIDIndex_ << " is out of bounds. The maximum value is " << parentageIDLookup_.size()-1 << ".\n"
1964  << "This should never happen.\n"
1965  << "Please report this to the framework hypernews forum 'hn-cms-edmFramework@cern.ch'.\n";
1966  }
1967  retValue.emplace(BranchID(prov.branchID_), parentageIDLookup_[prov.parentageIDIndex_]);
1968  }
1969  }
1970  return retValue;
1971  }
1972 
1974  public:
1975  explicit FullProvenanceReader(RootTree* rootTree, DaqProvenanceHelper const* daqProvenanceHelper);
1977  std::set<ProductProvenance> readProvenance(unsigned int transitionIndex) const override;
1978  private:
1979  void readProvenanceAsync(WaitingTask* task,
1980  ModuleCallingContext const* moduleCallingContext,
1981  unsigned int transitionIndex,
1982  std::atomic<const std::set<ProductProvenance>*>& writeTo
1983  ) const override;
1984 
1989  std::shared_ptr<std::recursive_mutex> mutex_;
1991  };
1992 
1995  rootTree_(rootTree),
1996  infoVector_(),
1997  pInfoVector_(&infoVector_),
1998  daqProvenanceHelper_(daqProvenanceHelper),
1999  mutex_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().second),
2000  acquirer_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().first)
2001 {
2002  }
2003 
2004  void
2006  ModuleCallingContext const* moduleCallingContext,
2007  unsigned int transitionIndex,
2008  std::atomic<const std::set<ProductProvenance>*>& writeTo
2009  ) const {
2010  readProvenanceAsyncImpl(this,
2012  task,
2013  transitionIndex,
2014  writeTo,
2015  moduleCallingContext,
2018  }
2019 
2020  std::set<ProductProvenance>
2021  FullProvenanceReader::readProvenance(unsigned int transitionIndex) const {
2022  {
2023  std::lock_guard<std::recursive_mutex> guard(*mutex_);
2025  }
2026  std::set<ProductProvenance> retValue;
2027  if(daqProvenanceHelper_) {
2028  for(auto const& info : infoVector_) {
2029  retValue.emplace(daqProvenanceHelper_->mapBranchID(info.branchID()),
2030  daqProvenanceHelper_->mapParentageID(info.parentageID()));
2031  }
2032  } else {
2033  for(auto const& info : infoVector_) {
2034  retValue.emplace(info);
2035  }
2036  }
2037  return retValue;
2038  }
2039 
2041  public:
2042  explicit OldProvenanceReader(RootTree* rootTree, EntryDescriptionMap const& theMap, DaqProvenanceHelper const* daqProvenanceHelper);
2043  ~OldProvenanceReader() override {}
2044  std::set<ProductProvenance> readProvenance(unsigned int transitionIndex) const override;
2045  private:
2046  void readProvenanceAsync(WaitingTask* task,
2047  ModuleCallingContext const* moduleCallingContext,
2048  unsigned int transitionIndex,
2049  std::atomic<const std::set<ProductProvenance>*>& writeTo
2050  ) const override;
2051 
2053  std::vector<EventEntryInfo> infoVector_;
2054  mutable std::vector<EventEntryInfo> *pInfoVector_;
2057  std::shared_ptr<std::recursive_mutex> mutex_;
2059  };
2060 
2061  OldProvenanceReader::OldProvenanceReader(RootTree* rootTree, EntryDescriptionMap const& theMap, DaqProvenanceHelper const* daqProvenanceHelper) :
2063  rootTree_(rootTree),
2064  infoVector_(),
2066  entryDescriptionMap_(theMap),
2067  daqProvenanceHelper_(daqProvenanceHelper),
2068  mutex_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().second),
2069  acquirer_(SharedResourcesRegistry::instance()->createAcquirerForSourceDelayedReader().first)
2070 {
2071  }
2072 
2073  void
2075  ModuleCallingContext const* moduleCallingContext,
2076  unsigned int transitionIndex,
2077  std::atomic<const std::set<ProductProvenance>*>& writeTo
2078  ) const {
2079  readProvenanceAsyncImpl(this,
2081  task,
2082  transitionIndex,
2083  writeTo,
2084  moduleCallingContext,
2085  rootTree_->rootDelayedReader()->preEventReadFromSourceSignal(),
2086  rootTree_->rootDelayedReader()->postEventReadFromSourceSignal());
2087  }
2088 
2089  std::set<ProductProvenance>
2090  OldProvenanceReader::readProvenance(unsigned int transitionIndex) const {
2091  {
2092  std::lock_guard<std::recursive_mutex> guard(*mutex_);
2093  rootTree_->branchEntryInfoBranch()->SetAddress(&pInfoVector_);
2094  roottree::getEntry(rootTree_->branchEntryInfoBranch(), rootTree_->entryNumberForIndex(transitionIndex));
2095  }
2096  std::set<ProductProvenance> retValue;
2097  for(auto const& info : infoVector_) {
2098  EntryDescriptionMap::const_iterator iter = entryDescriptionMap_.find(info.entryDescriptionID());
2099  assert(iter != entryDescriptionMap_.end());
2100  Parentage parentage(iter->second.parents());
2101  if(daqProvenanceHelper_) {
2102  retValue.emplace(daqProvenanceHelper_->mapBranchID(info.branchID()),
2103  daqProvenanceHelper_->mapParentageID(parentage.id()));
2104  } else {
2105  retValue.emplace(info.branchID(), parentage.id());
2106  }
2107  }
2108  return retValue;
2109  }
2110 
2112  public:
2115  private:
2116  std::set<ProductProvenance> readProvenance(unsigned int) const override;
2117  void readProvenanceAsync(WaitingTask* task,
2118  ModuleCallingContext const* moduleCallingContext,
2119  unsigned int transitionIndex,
2120  std::atomic<const std::set<ProductProvenance>*>& writeTo
2121  ) const override;
2122  };
2123 
2126  }
2127 
2128  std::set<ProductProvenance>
2130  // Not providing parentage!!!
2131  return std::set<ProductProvenance>{};
2132  }
2133  void
2135  ModuleCallingContext const* moduleCallingContext,
2136  unsigned int transitionIndex,
2137  std::atomic<const std::set<ProductProvenance>*>& writeTo
2138  ) const {
2139  if(nullptr == writeTo.load()) {
2140  auto emptyProv = std::make_unique<const std::set<ProductProvenance>>();
2141  const std::set<ProductProvenance>* expected = nullptr;
2142  if(writeTo.compare_exchange_strong(expected,emptyProv.get())) {
2143  emptyProv.release();
2144  }
2145  }
2146  }
2147 
2148  std::unique_ptr<ProvenanceReaderBase>
2150  return std::make_unique<DummyProvenanceReader>();
2151  }
2152 
2153  std::unique_ptr<ProvenanceReaderBase>
2154  MakeOldProvenanceReader::makeReader(RootTree& rootTree, DaqProvenanceHelper const* daqProvenanceHelper) const {
2155  return std::make_unique<OldProvenanceReader>(&rootTree, *entryDescriptionMap_, daqProvenanceHelper);
2156  }
2157 
2158  std::unique_ptr<ProvenanceReaderBase>
2159  MakeFullProvenanceReader::makeReader(RootTree& rootTree, DaqProvenanceHelper const* daqProvenanceHelper) const {
2160  return std::make_unique<FullProvenanceReader>(&rootTree, daqProvenanceHelper);
2161  }
2162 
2163  std::unique_ptr<ProvenanceReaderBase>
2164  MakeReducedProvenanceReader::makeReader(RootTree& rootTree, DaqProvenanceHelper const* daqProvenanceHelper) const {
2165  return std::make_unique<ReducedProvenanceReader>(&rootTree, parentageIDLookup_, daqProvenanceHelper);
2166  }
2167 }
void dropOnInput(ProductRegistry &reg, ProductSelectorRules const &rules, bool dropDescendants, InputType inputType)
Definition: RootFile.cc:1685
RunNumber_t run() const
Definition: EventID.h:39
EventID const & eventID() const
Definition: RootFile.h:172
Int_t getEntry(TBranch *branch, EntryNumber entryNumber)
Definition: RootTree.cc:510
void fillEventNumbersOrEntries(bool needEventNumbers, bool needEventEntries) const
EventNumber_t event() const
Definition: EventID.h:41
std::string const & idToParameterSetBlobsBranchName()
Definition: BranchType.cc:255
edm::propagate_const< std::unique_ptr< ThinnedAssociationsHelper > > fileThinnedAssociationsHelper_
Definition: RootFile.h:284
std::vector< ProcessConfiguration > ProcessConfigurationVector
void setEventFinder(std::shared_ptr< EventFinder > ptr) const
std::string const & branchName() const
StoredProductProvenanceVector provVector_
Definition: RootFile.cc:1859
bool isRealData() const
bool selected(BranchDescription const &desc) const
std::string const & BranchTypeToAuxiliaryBranchName(BranchType const &branchType)
Definition: BranchType.cc:115
unsigned int const defaultNonEventLearningEntries
Definition: RootTree.h:49
InputType
Definition: InputType.h:5
static const TGPicture * info(bool iBackgroundIsBlack)
void doneFileInitialization() const
Clears the temporary vector of event numbers to reduce memory usage.
static constexpr int invalidIndex
TPRegexp parents
Definition: eve_filter.cc:21
BranchType const & branchType() const
std::string const & parentageTreeName()
Definition: BranchType.cc:159
std::vector< BranchIDList > BranchIDLists
Definition: BranchIDList.h:19
std::shared_ptr< BranchIDLists const > branchIDLists_
Definition: RootFile.h:282
std::string const & entryDescriptionBranchName()
Definition: BranchType.cc:154
bool skipEvents(int &offset)
Definition: RootFile.cc:1308
StreamContext const * getStreamContext() const
edm::propagate_const< std::shared_ptr< DuplicateChecker > > duplicateChecker_
Definition: RootFile.h:294
FileFormatVersion fileFormatVersion() const
Definition: RootFile.h:176
std::shared_ptr< LuminosityBlockAuxiliary > fillLumiAuxiliary()
Definition: RootFile.cc:1263
ProductProvenanceVector infoVector_
Definition: RootFile.cc:1986
static std::string const source("source")
static Timestamp invalidTimestamp()
Definition: Timestamp.h:101
ForwardSequence::const_iterator lower_bound_all(ForwardSequence const &s, Datum const &d)
wrappers for std::lower_bound
Definition: Algorithms.h:91
Definition: chain.py:1
void readEntryDescriptionTree(EntryDescriptionMap &entryDescriptionMap, InputType inputType)
Definition: RootFile.cc:533
bool branchListIndexesUnchanged() const
Definition: RootFile.h:179
RunNumber_t run() const
Definition: RunID.h:39
std::vector< std::string > const & branchNames() const
Definition: RootTree.h:136
bool setEntryAtNextEventInLumi(RunNumber_t run, LuminosityBlockNumber_t lumi)
Definition: RootFile.cc:1633
RootFile(std::string const &fileName, ProcessConfiguration const &processConfiguration, std::string const &logicalFileName, std::shared_ptr< InputFile > filePtr, std::shared_ptr< EventSkipperByID > eventSkipperByID, bool skipAnyEvents, int remainingEvents, int remainingLumis, unsigned int nStreams, unsigned int treeCacheSize, int treeMaxVirtualSize, InputSource::ProcessingMode processingMode, RunHelperBase *runHelper, bool noEventSort, ProductSelectorRules const &productSelectorRules, InputType inputType, std::shared_ptr< BranchIDListHelper > branchIDListHelper, std::shared_ptr< ThinnedAssociationsHelper > thinnedAssociationsHelper, std::vector< BranchID > const *associationsFromSecondary, std::shared_ptr< DuplicateChecker > duplicateChecker, bool dropDescendantsOfDroppedProducts, ProcessHistoryRegistry &processHistoryRegistry, std::vector< std::shared_ptr< IndexIntoFile > > const &indexesIntoFiles, std::vector< std::shared_ptr< IndexIntoFile > >::size_type currentIndexIntoFile, std::vector< ProcessHistoryID > &orderedProcessHistoryIDs, bool bypassVersionCheck, bool labelRawDataLikeMC, bool usingGoToEvent, bool enablePrefetching)
Definition: RootFile.cc:129
static PFTauRenderPlugin instance
MakeOldProvenanceReader(std::unique_ptr< EntryDescriptionMap > &&entryDescriptionMap)
Definition: RootFile.cc:73
int whyNotFastClonable_
Definition: RootFile.h:272
std::vector< BranchID > & parentsForUpdate()
Definition: Parentage.h:45
edm::propagate_const< std::shared_ptr< EventSkipperByID > > eventSkipperByID_
Definition: RootFile.h:258
edm::propagate_const< std::shared_ptr< ThinnedAssociationsHelper > > thinnedAssociationsHelper_
Definition: RootFile.h:285
ParentageID id() const
Definition: Parentage.cc:23
DaqProvenanceHelper const * daqProvenanceHelper_
Definition: RootFile.cc:1988
std::map< std::string, std::string > newBranchToOldBranch_
Definition: RootFile.h:288
bool empty() const
True if no runs, lumis, or events are in the file.
std::map< BranchKey, BranchDescription > ProductList
BranchID const & mapBranchID(BranchID const &branchID) const
bool registerProcessHistory(ProcessHistory const &processHistory)
bool setEntryAtRun(RunNumber_t run)
Definition: RootFile.cc:1625
void readRun_(RunPrincipal &runPrincipal)
Definition: RootFile.cc:1532
ProductProvenanceVector * pInfoVector_
Definition: RootFile.cc:1987
RootTree lumiTree_
Definition: RootFile.h:277
Timestamp const & time() const
edm::propagate_const< ProcessHistoryRegistry * > processHistoryRegistry_
Definition: RootFile.h:256
edm::propagate_const< std::unique_ptr< DaqProvenanceHelper > > daqProvenanceHelper_
Definition: RootFile.h:299
edm::propagate_const< std::shared_ptr< RunAuxiliary > > savedRunAuxiliary_
Definition: RootFile.h:269
unsigned long long EventNumber_t
EntryNumber const & entries() const
Definition: RootTree.h:133
RunNumber_t run() const
void setSignals(signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const *preEventReadSource, signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const *postEventReadSource)
Definition: RootTree.cc:501
std::string const & fileFormatVersionBranchName()
Definition: BranchType.cc:218
IndexIntoFileItr begin(SortOrder sortOrder) const
FullProvenanceReader(RootTree *rootTree, DaqProvenanceHelper const *daqProvenanceHelper)
Definition: RootFile.cc:1993
EntryDescriptionMap const & entryDescriptionMap_
Definition: RootFile.cc:2055
std::vector< BranchID > const & parents() const
std::string const & eventSelectionsBranchName()
Definition: BranchType.cc:243
edm::propagate_const< std::unique_ptr< EntryDescriptionMap > > entryDescriptionMap_
Definition: RootFile.cc:76
void push(T &&iAction)
asynchronously pushes functor iAction into queue
InputSource::ProcessingMode processingMode_
Definition: RootFile.h:286
~RootFileEventFinder() override
Definition: RootFile.cc:112
IndexIntoFile::IndexIntoFileItr indexIntoFileIter_
Definition: RootFile.h:266
std::shared_ptr< ProductProvenanceRetriever const > eventProductProvenanceRetriever(size_t index) const
Definition: RootFile.h:250
bool current() const
Definition: RootTree.h:126
void fillRunPrincipal(ProcessHistoryRegistry const &processHistoryRegistry, DelayedReader *reader=0)
Definition: RunPrincipal.cc:20
void stable_sort_all(RandomAccessSequence &s)
wrappers for std::stable_sort
Definition: Algorithms.h:135
std::shared_ptr< std::recursive_mutex > mutex_
Definition: RootFile.cc:2057
#define nullptr
LuminosityBlockNumber_t lumi() const
void fillBranchEntryMeta(TBranch *branch, T *&pbuf)
Definition: RootTree.h:145
std::unique_ptr< FileBlock > createFileBlock() const
Definition: RootFile.cc:661
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
void trainCache(char const *branchNames)
Definition: RootTree.cc:471
bool readCurrentEvent(EventPrincipal &cache)
Definition: RootFile.cc:1433
SharedResourcesAcquirer acquirer_
Definition: RootFile.cc:2058
~FullProvenanceReader() override
Definition: RootFile.cc:1976
ServiceToken presentToken() const
std::shared_ptr< ProductRegistry const > productRegistry() const
Definition: RootFile.h:167
uint16_t size_type
unsigned int LuminosityBlockNumber_t
IndexIntoFileItr findRunPosition(RunNumber_t run) const
Same as findPosition.
void readProvenanceAsync(WaitingTask *task, ModuleCallingContext const *moduleCallingContext, unsigned int transitionIndex, std::atomic< const std::set< ProductProvenance > * > &writeTo) const override
Definition: RootFile.cc:2074
TBranch * branchEntryInfoBranch() const
Definition: RootTree.h:184
bool fillEventAuxiliary(IndexIntoFile::EntryNumber_t entry)
Definition: RootFile.cc:1195
MakeReducedProvenanceReader(std::vector< ParentageID > const &parentageIDLookup)
Definition: RootFile.cc:84
IndexIntoFileItr findPosition(RunNumber_t run, LuminosityBlockNumber_t lumi=0U, EventNumber_t event=0U) const
std::vector< EventSelectionID > EventSelectionIDVector
edm::propagate_const< TClass * > edProductClass_
Definition: RootFile.h:300
LuminosityBlockNumber_t luminosityBlock() const
FileFormatVersion fileFormatVersion_
Definition: RootFile.h:259
TTree const * metaTree() const
Definition: RootTree.h:180
std::string const & parameterSetsTreeName()
Definition: BranchType.cc:251
void setPosition(IndexIntoFile::IndexIntoFileItr const &position)
Definition: RootFile.cc:692
void readAllFromSourceAndMergeImmediately()
Definition: Principal.cc:916
std::vector< edm::propagate_const< std::shared_ptr< ProductProvenanceRetriever > > > eventProductProvenanceRetrievers_
Definition: RootFile.h:297
edm::propagate_const< RootTree * > rootTree_
Definition: RootFile.cc:2052
ProductList const & productList() const
U second(std::pair< T, U > const &p)
std::vector< EventProcessHistoryID >::const_iterator eventProcessHistoryIter_
Definition: RootFile.h:268
void reduceProcessHistoryIDs(ProcessHistoryRegistry const &processHistoryRegistry)
bool skipThisEntry()
Definition: RootFile.cc:702
TTree const * tree() const
Definition: RootTree.h:178
FileID fid_
Definition: RootFile.h:260
long long EntryNumber_t
bool containsItem(RunNumber_t run, LuminosityBlockNumber_t lumi, EventNumber_t event) const
Definition: RootFile.cc:757
void setParents(std::vector< BranchID > const &parents)
Definition: Parentage.h:46
DaqProvenanceHelper const * daqProvenanceHelper_
Definition: RootFile.cc:1862
std::string const logicalFile_
Definition: RootFile.h:254
def principal(options)
std::vector< BranchListIndex > BranchListIndexes
IndexIntoFile::IndexIntoFileItr indexIntoFileEnd_
Definition: RootFile.h:265
std::string moduleName(Provenance const &provenance)
Definition: Provenance.cc:27
void fillLuminosityBlockPrincipal(ProcessHistoryRegistry const &processHistoryRegistry, DelayedReader *reader=0)
void readParentageTree(InputType inputType)
Definition: RootFile.cc:579
std::string const & processHistoryMapBranchName()
Definition: BranchType.cc:193
std::vector< ProcessHistoryID > & setProcessHistoryIDs()
bool next()
Definition: RootTree.h:123
bool noEventSort_
Definition: RootFile.h:271
std::string const & className() const
bool wasLastEventJustRead() const
Definition: RootFile.cc:808
void insertEntryForIndex(unsigned int index)
Definition: RootTree.cc:100
DelayedReader * resetAndGetRootDelayedReader() const
Definition: RootTree.cc:118
bool setEntryAtLumi(RunNumber_t run, LuminosityBlockNumber_t lumi)
Definition: RootFile.cc:1617
std::string friendlyName(std::string const &iFullName)
static constexpr RunNumber_t invalidRun
std::string const & entryDescriptionTreeName()
Definition: BranchType.cc:146
void close()
Definition: RootFile.cc:1164
virtual signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const * postEventReadFromSourceSignal() const =0
std::string const & fid() const
Definition: FileID.h:19
std::shared_ptr< std::recursive_mutex > mutex_
Definition: RootFile.cc:1989
RootTree const & runTree() const
Definition: RootFile.h:175
bool modifiedIDs() const
Definition: RootFile.h:180
IndexIntoFile::EntryNumber_t lastEventEntryNumberRead_
Definition: RootFile.h:280
void copyPosition(IndexIntoFileItr const &position)
Copy the position without modifying the pointer to the IndexIntoFile or size.
std::string const & indexIntoFileBranchName()
Definition: BranchType.cc:233
IndexIntoFileItr end(SortOrder sortOrder) const
Used to end an iteration over the Runs, Lumis, and Events in a file.
std::string const & eventHistoryBranchName()
Definition: BranchType.cc:238
std::set< ProductProvenance > readProvenance(unsigned int transitionIndex) const override
Definition: RootFile.cc:2021
void skipEventBackward(int &phIndexOfEvent, RunNumber_t &runOfEvent, LuminosityBlockNumber_t &lumiOfEvent, EntryNumber_t &eventEntry)
Long64_t numEntries(TFile *hdl, std::string const &trname)
Definition: CollUtil.cc:50
edm::propagate_const< RootTree * > rootTree_
Definition: RootFile.cc:1857
std::shared_ptr< RunAuxiliary > readRunAuxiliary_()
Definition: RootFile.cc:1466
EventAuxiliary const & eventAux() const
Definition: RootFile.h:168
EventSelectionIDVector eventSelectionIDs_
Definition: RootFile.h:290
std::shared_ptr< LuminosityBlockAuxiliary > readLuminosityBlockAuxiliary_()
Definition: RootFile.cc:1552
bool eventHistoryTree() const
bool isEarlierRelease(std::string const &a, std::string const &b)
StoredProductProvenanceVector const * pProvVector_
Definition: RootFile.cc:1860
int whyNotFastClonable() const
Definition: RootFile.h:177
edm::propagate_const< std::unique_ptr< ProvenanceAdaptor > > provenanceAdaptor_
Definition: RootFile.h:295
StreamID streamID() const
static ServiceRegistry & instance()
RunNumber_t run() const
std::shared_ptr< RunAuxiliary > fillRunAuxiliary()
Definition: RootFile.cc:1287
bool skipAnyEvents_
Definition: RootFile.h:270
void setSignals(signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const *preEventReadSource, signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const *postEventReadSource)
Definition: RootFile.cc:1809
RootTree eventTree_
Definition: RootFile.h:276
std::vector< BranchDescription const * > allBranchDescriptions() const
EventNumber_t getEventNumberOfEntry(roottree::EntryNumber entry) const override
Definition: RootFile.cc:114
std::map< ParameterSetID, ParameterSetID > ParameterSetIdConverter
bool isValid() const
Definition: FileID.h:18
bool iterationWillBeInEntryOrder(SortOrder sortOrder) const
Used to determine whether or not to disable fast cloning.
std::string const & friendlyClassName() const
BranchID const & branchID() const
TypeWithDict const & unwrappedType() const
bool insertMapped(value_type const &v, bool forceUpdate=false)
Definition: Registry.cc:36
std::array< bool, NumBranchTypes > const & hasNewlyDroppedBranch() const
Definition: RootFile.h:178
IndexIntoFile::EntryType getNextItemType(RunNumber_t &run, LuminosityBlockNumber_t &lumi, EventNumber_t &event)
Definition: RootFile.cc:762
BranchListIndexes branchListIndexes_
Definition: RootFile.h:291
std::string const & metaDataTreeName()
Definition: BranchType.cc:168
virtual signalslot::Signal< void(StreamContext const &, ModuleCallingContext const &)> const * preEventReadFromSourceSignal() const =0
bool isDuplicateEvent()
Definition: RootFile.cc:746
LuminosityBlockNumber_t oldLuminosityBlock() const
EntryDescriptionID id() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
IndexIntoFile::EntryNumber_t EntryNumber
Definition: RootTree.h:50
bool wasFirstEventJustRead() const
Definition: RootFile.cc:815
std::string const & parameterSetMapBranchName()
Definition: BranchType.cc:183
Hash< ProcessHistoryType > ProcessHistoryID
std::unique_ptr< ProvenanceReaderBase > makeReader(RootTree &eventTree, DaqProvenanceHelper const *daqProvenanceHelper) const override
Definition: RootFile.cc:2149
edm::propagate_const< std::unique_ptr< History > > history_
Definition: RootFile.h:292
SerialTaskQueueChain & serialQueueChain() const
std::vector< ProductProvenance > ProductProvenanceVector
const int drop
std::string const & processHistoryBranchName()
Definition: BranchType.cc:198
IndexIntoFile & indexIntoFile_
Definition: RootFile.h:262
void readEvent(EventPrincipal &cache)
Definition: RootFile.cc:1410
RootTree runTree_
Definition: RootFile.h:278
SharedResourcesAcquirer acquirer_
Definition: RootFile.cc:1990
edm::propagate_const< std::unique_ptr< MakeProvenanceReader > > provenanceReaderMaker_
Definition: RootFile.h:296
EntryNumber const & entryNumber() const
Definition: RootTree.h:131
bool containsItem(RunNumber_t run, LuminosityBlockNumber_t lumi, EventNumber_t event) const
bool goToEvent(EventID const &eventID)
Definition: RootFile.cc:1376
std::string getReleaseVersion()
std::shared_ptr< RunAuxiliary const > savedRunAuxiliary() const
Definition: RootFile.h:244
std::set< ProductProvenance > readProvenance(unsigned int) const override
Definition: RootFile.cc:1945
void setAtEventEntry(IndexIntoFile::EntryNumber_t entry)
Definition: RootFile.cc:1461
bool storedProductProvenanceUsed() const
void fixIndexes(std::vector< ProcessHistoryID > &processHistoryIDs)
DelayedReader * rootDelayedReader() const
Definition: RootTree.cc:124
std::vector< RunOrLumiEntry > & setRunOrLumiEntries()
EntryNumber const & entryNumberForIndex(unsigned int index) const
Definition: RootTree.cc:94
void conversion(EventAux const &from, EventAuxiliary &to)
Definition: EventAux.cc:9
std::array< bool, NumBranchTypes > hasNewlyDroppedBranch_
Definition: RootFile.h:273
std::vector< EventNumber_t > & unsortedEventNumbers() const
unsigned int value() const
Definition: StreamID.h:46
static constexpr EntryNumber_t invalidEntry
ForwardSequence::const_iterator find_in_all(ForwardSequence const &s, Datum const &d)
wrappers for std::find
Definition: Algorithms.h:32
std::shared_ptr< ProductRegistry const > productRegistry_
Definition: RootFile.h:281
void setIfFastClonable(int remainingEvents, int remainingLumis)
Definition: RootFile.cc:615
std::unique_ptr< ProvenanceReaderBase > makeReader(RootTree &eventTree, DaqProvenanceHelper const *daqProvenanceHelper) const override
Definition: RootFile.cc:2164
void fillBranchEntry(TBranch *branch, T *&pbuf)
Definition: RootTree.h:156
static constexpr LuminosityBlockNumber_t invalidLumi
IndexIntoFileItr findLumiPosition(RunNumber_t run, LuminosityBlockNumber_t lumi) const
std::map< EntryDescriptionID, EventEntryDescription > EntryDescriptionMap
Definition: RootFile.h:51
edm::propagate_const< std::shared_ptr< BranchChildren > > branchChildren_
Definition: RootFile.h:293
std::string const & parentageBranchName()
Definition: BranchType.cc:163
edm::propagate_const< std::shared_ptr< InputFile > > filePtr_
Definition: RootFile.h:257
std::shared_ptr< BranchChildren const > branchChildren() const
Definition: RootFile.h:247
double b
Definition: hdecay.h:120
LuminosityBlockNumber_t luminosityBlock() const
bool processHistorySameWithinRun() const
void fillEventHistory()
Definition: RootFile.cc:1205
ProcessHistoryID const & processHistoryID(int i) const
ProductList & productListUpdator()
void setProcessHistoryID(ProcessHistoryID const &phid)
unsigned int transitionIndex() const
void readProvenanceAsync(WaitingTask *task, ModuleCallingContext const *moduleCallingContext, unsigned int transitionIndex, std::atomic< const std::set< ProductProvenance > * > &writeTo) const override
Definition: RootFile.cc:2134
void readProvenanceAsync(WaitingTask *task, ModuleCallingContext const *moduleCallingContext, unsigned int transitionIndex, std::atomic< const std::set< ProductProvenance > * > &writeTo) const override
Definition: RootFile.cc:1929
std::vector< ParentageID > const & parentageIDLookup_
Definition: RootFile.cc:1861
std::shared_ptr< std::recursive_mutex > mutex_
Definition: RootFile.cc:1863
std::unique_ptr< ProvenanceReaderBase > makeReader(RootTree &eventTree, DaqProvenanceHelper const *daqProvenanceHelper) const override
Definition: RootFile.cc:2154
bool isValid() const
Definition: RootTree.cc:106
~OldProvenanceReader() override
Definition: RootFile.cc:2043
std::shared_ptr< ProductProvenanceRetriever > makeProductProvenanceRetriever(unsigned int iStreamIndex)
Definition: RootFile.cc:1833
std::vector< StoredProductProvenance > StoredProductProvenanceVector
std::string const & file() const
Definition: RootFile.h:166
std::string const & productDescriptionBranchName()
Definition: BranchType.cc:173
std::string const & processConfigurationBranchName()
Definition: BranchType.cc:203
IndexIntoFile::IndexIntoFileItr indexIntoFileBegin_
Definition: RootFile.h:264
ProcessHistoryID const & processHistoryID() const
EventID const & id() const
edm::propagate_const< RunHelperBase * > runHelper_
Definition: RootFile.h:287
bool branchListIndexesUnchanged_
Definition: RootFile.h:274
void fillEventPrincipal(EventAuxiliary const &aux, ProcessHistoryRegistry const &processHistoryRegistry, DelayedReader *reader=0)
void fillIndexIntoFile()
Definition: RootFile.cc:877
#define begin
Definition: vmac.h:32
HLT enums.
edm::propagate_const< TTree * > eventHistoryTree_
Definition: RootFile.h:289
void readLuminosityBlock_(LuminosityBlockPrincipal &lumiPrincipal)
Definition: RootFile.cc:1590
std::vector< EventEntryInfo > * pInfoVector_
Definition: RootFile.cc:2054
std::unique_ptr< ProvenanceReaderBase > makeReader(RootTree &eventTree, DaqProvenanceHelper const *daqProvenanceHelper) const override
Definition: RootFile.cc:2159
unsigned int const defaultNonEventCacheSize
Definition: RootTree.h:47
RootTreePtrArray treePointers_
Definition: RootFile.h:279
void fillEventNumbers() const
double a
Definition: hdecay.h:121
void initialize(ProductSelectorRules const &rules, std::vector< BranchDescription const * > const &branchDescriptions)
static int position[264][3]
Definition: ReadPGInfo.cc:509
void initializeDuplicateChecker(std::vector< std::shared_ptr< IndexIntoFile > > const &indexesIntoFiles, std::vector< std::shared_ptr< IndexIntoFile > >::size_type currentIndexIntoFile)
Definition: RootFile.cc:1659
bool perEventProductIDs() const
std::string const & BranchTypeToProductProvenanceBranchName(BranchType const &BranchType)
Definition: BranchType.cc:131
std::vector< ProcessHistoryID > & orderedProcessHistoryIDs_
Definition: RootFile.h:263
std::string const & productDependenciesBranchName()
Definition: BranchType.cc:178
void setNumberOfEvents(EntryNumber_t nevents) const
void initAssociationsFromSecondary(std::vector< BranchID > const &)
Definition: RootFile.cc:697
std::string const & thinnedAssociationsHelperBranchName()
Definition: BranchType.cc:213
std::string const & entryDescriptionIDBranchName()
Definition: BranchType.cc:150
~DummyProvenanceReader() override
Definition: RootFile.cc:2114
bool useReducedProcessHistoryID() const
edm::propagate_const< TBranch * > provBranch_
Definition: RootFile.cc:1858
std::unique_ptr< MakeProvenanceReader > makeProvenanceReaderMaker(InputType inputType)
Definition: RootFile.cc:1816
std::unique_ptr< WrapperBase > getWrapperBasePtr(void *p, int offset)
std::string const & branchIDListBranchName()
Definition: BranchType.cc:208
BranchID const & originalBranchID() const
unsigned int RunNumber_t
std::string const & branchListIndexesBranchName()
Definition: BranchType.cc:247
std::vector< ParentageID > parentageIDLookup_
Definition: RootFile.h:298
EventAuxiliary eventAux_
Definition: RootFile.h:275
void readProvenanceAsync(WaitingTask *task, ModuleCallingContext const *moduleCallingContext, unsigned int transitionIndex, std::atomic< const std::set< ProductProvenance > * > &writeTo) const override
Definition: RootFile.cc:2005
void fillThisEventAuxiliary()
Definition: RootFile.cc:1176
void resetTraining()
Definition: RootTree.h:190
bool setEntryAtEvent(RunNumber_t run, LuminosityBlockNumber_t lumi, EventNumber_t event)
Definition: RootFile.cc:1609
std::vector< ParentageID > const & parentageIDLookup_
Definition: RootFile.cc:87
ReducedProvenanceReader(RootTree *iRootTree, std::vector< ParentageID > const &iParentageIDLookup, DaqProvenanceHelper const *daqProvenanceHelper)
Definition: RootFile.cc:1867
std::string const & eventHistoryTreeName()
Definition: BranchType.cc:268
static Interceptor::Registry registry("Interceptor")
RootFileEventFinder(RootTree &eventTree)
Definition: RootFile.cc:111
ParentageID const & mapParentageID(ParentageID const &phid) const
edm::propagate_const< std::shared_ptr< BranchIDListHelper > > branchIDListHelper_
Definition: RootFile.h:283
void validateFile(InputType inputType, bool usingGoToEvent)
Definition: RootFile.cc:1103
IndexIntoFile::IndexIntoFileItr indexIntoFileIter() const
Definition: RootFile.cc:687
std::string const & newBranchToOldBranch(std::string const &newBranch) const
Definition: RootFile.cc:678
T first(std::pair< T, U > const &p)
std::unique_ptr< TTreeCache > trainCache(TTree *tree, InputFile &file, unsigned int cacheSize, char const *branchNames)
Definition: RootTree.cc:534
static ParentageRegistry * instance()
OldProvenanceReader(RootTree *rootTree, EntryDescriptionMap const &theMap, DaqProvenanceHelper const *daqProvenanceHelper)
Definition: RootFile.cc:2061
SharedResourcesAcquirer acquirer_
Definition: RootFile.cc:1864
IndexIntoFileItr findEventPosition(RunNumber_t run, LuminosityBlockNumber_t lumi, EventNumber_t event) const
unsigned int const defaultLearningEntries
Definition: RootTree.h:48
std::vector< EventEntryInfo > infoVector_
Definition: RootFile.cc:2053
void markBranchToBeDropped(bool dropDescendants, BranchDescription const &branch, std::set< BranchID > &branchesToDrop, std::map< BranchID, BranchID > const &droppedToKeptAlias) const
Definition: RootFile.cc:1675
LuminosityBlockNumber_t peekAheadAtLumi() const
bool hasThinnedAssociations() const
std::string const & fileIdentifierBranchName()
Definition: BranchType.cc:223
std::string const & wrappedName() const
EventNumber_t event() const
std::string const & moduleDescriptionMapBranchName()
Definition: BranchType.cc:188
bool insertMapped(value_type const &v)
std::string const file_
Definition: RootFile.h:253
def move(src, dest)
Definition: eostools.py:510
void skipEventForward(int &phIndexOfSkippedEvent, RunNumber_t &runOfSkippedEvent, LuminosityBlockNumber_t &lumiOfSkippedEvent, EntryNumber_t &skippedEventEntry)
std::vector< EventProcessHistoryID > eventProcessHistoryIDs_
Definition: RootFile.h:267
static Registry * instance()
Definition: Registry.cc:13
std::set< ProductProvenance > readProvenance(unsigned int transitionIndex) const override
Definition: RootFile.cc:2090
DaqProvenanceHelper const * daqProvenanceHelper_
Definition: RootFile.cc:2056
std::string createGlobalIdentifier()
Definition: event.py:1
void copyProduct(BranchDescription const &productdesc)
void readEventHistoryTree()
Definition: RootFile.cc:1646
void fillAux(T *&pAux)
Definition: RootTree.h:140
std::set< ProductProvenance > readProvenance(unsigned int) const override
Definition: RootFile.cc:2129
void setEntryNumber(EntryNumber theEntryNumber)
Definition: RootTree.cc:205
void reportOpened(std::string const &inputType)
Definition: RootFile.cc:1150
def operate(timelog, memlog, json_f, num)