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