Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-15 02:27:58

0001 // -*- C++ -*-
0002 //
0003 // Package:     Services
0004 // Class  :     Tracer
0005 //
0006 // Implementation:
0007 //     <Notes on implementation>
0008 //
0009 // Original Author:  Chris Jones
0010 //         Created:  Thu Sep  8 14:17:58 EDT 2005
0011 //
0012 
0013 #include "FWCore/Framework/interface/ComponentDescription.h"
0014 #include "FWCore/Framework/interface/DataKey.h"
0015 #include "FWCore/Framework/interface/ESProducer.h"
0016 #include "FWCore/Framework/interface/ESProductResolverProvider.h"
0017 #include "FWCore/Framework/interface/EventSetupRecordKey.h"
0018 #include "FWCore/Framework/interface/IOVSyncValue.h"
0019 
0020 #include "FWCore/ServiceRegistry/interface/ServiceMaker.h"
0021 
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/ServiceRegistry/interface/ActivityRegistry.h"
0024 
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 #include "FWCore/ServiceRegistry/interface/ModuleConsumesESInfo.h"
0027 #include "FWCore/ServiceRegistry/interface/ModuleConsumesInfo.h"
0028 #include "FWCore/ServiceRegistry/interface/ESModuleConsumesInfo.h"
0029 #include "FWCore/ServiceRegistry/interface/PathsAndConsumesOfModulesBase.h"
0030 #include "FWCore/ServiceRegistry/interface/SystemBounds.h"
0031 #include "FWCore/Utilities/interface/BranchType.h"
0032 #include "FWCore/Utilities/interface/Exception.h"
0033 #include "FWCore/Utilities/interface/ProductKindOfType.h"
0034 #include "FWCore/Utilities/interface/TimeOfDay.h"
0035 
0036 #include "DataFormats/Provenance/interface/EventID.h"
0037 #include "DataFormats/Provenance/interface/LuminosityBlockID.h"
0038 #include "DataFormats/Provenance/interface/ModuleDescription.h"
0039 #include "DataFormats/Provenance/interface/RunID.h"
0040 #include "DataFormats/Provenance/interface/Timestamp.h"
0041 
0042 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0043 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0044 #include "FWCore/ServiceRegistry/interface/GlobalContext.h"
0045 #include "FWCore/ServiceRegistry/interface/ModuleCallingContext.h"
0046 #include "FWCore/ServiceRegistry/interface/ESModuleCallingContext.h"
0047 #include "FWCore/ServiceRegistry/interface/PathContext.h"
0048 #include "FWCore/ServiceRegistry/interface/ProcessContext.h"
0049 #include "FWCore/ServiceRegistry/interface/StreamContext.h"
0050 #include "DataFormats/Common/interface/HLTPathStatus.h"
0051 
0052 #include "tracer_setupFile.h"
0053 
0054 #include <iostream>
0055 #include <vector>
0056 #include <string>
0057 #include <set>
0058 #include <optional>
0059 
0060 namespace edm {
0061   namespace service {
0062     class Tracer {
0063     public:
0064       Tracer(const ParameterSet&, ActivityRegistry&);
0065 
0066       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0067 
0068       void preallocate(service::SystemBounds const&);
0069 
0070       void preBeginJob(ProcessContext const&);
0071       void postBeginJob();
0072       void preEndJob();
0073       void postEndJob();
0074       void lookupInitializationComplete(PathsAndConsumesOfModulesBase const&, ProcessContext const&);
0075 
0076       void preBeginStream(StreamContext const&);
0077       void postBeginStream(StreamContext const&);
0078       void preEndStream(StreamContext const&);
0079       void postEndStream(StreamContext const&);
0080 
0081       void preSourceEvent(StreamID);
0082       void postSourceEvent(StreamID);
0083 
0084       void preSourceLumi(LuminosityBlockIndex);
0085       void postSourceLumi(LuminosityBlockIndex);
0086 
0087       void preSourceRun(RunIndex);
0088       void postSourceRun(RunIndex);
0089 
0090       void preSourceProcessBlock();
0091       void postSourceProcessBlock(std::string const&);
0092 
0093       void preOpenFile(std::string const&);
0094       void postOpenFile(std::string const&);
0095 
0096       void preCloseFile(std::string const& lfn);
0097       void postCloseFile(std::string const&);
0098 
0099       void preModuleBeginStream(StreamContext const&, ModuleCallingContext const&);
0100       void postModuleBeginStream(StreamContext const&, ModuleCallingContext const&);
0101 
0102       void preModuleEndStream(StreamContext const&, ModuleCallingContext const&);
0103       void postModuleEndStream(StreamContext const&, ModuleCallingContext const&);
0104 
0105       void preBeginProcessBlock(GlobalContext const&);
0106       void postBeginProcessBlock(GlobalContext const&);
0107 
0108       void preAccessInputProcessBlock(GlobalContext const&);
0109       void postAccessInputProcessBlock(GlobalContext const&);
0110 
0111       void preEndProcessBlock(GlobalContext const&);
0112       void postEndProcessBlock(GlobalContext const&);
0113 
0114       void preWriteProcessBlock(GlobalContext const&);
0115       void postWriteProcessBlock(GlobalContext const&);
0116 
0117       void preGlobalBeginRun(GlobalContext const&);
0118       void postGlobalBeginRun(GlobalContext const&);
0119 
0120       void preGlobalEndRun(GlobalContext const&);
0121       void postGlobalEndRun(GlobalContext const&);
0122 
0123       void preGlobalWriteRun(GlobalContext const&);
0124       void postGlobalWriteRun(GlobalContext const&);
0125 
0126       void preStreamBeginRun(StreamContext const&);
0127       void postStreamBeginRun(StreamContext const&);
0128 
0129       void preStreamEndRun(StreamContext const&);
0130       void postStreamEndRun(StreamContext const&);
0131 
0132       void preGlobalBeginLumi(GlobalContext const&);
0133       void postGlobalBeginLumi(GlobalContext const&);
0134 
0135       void preGlobalEndLumi(GlobalContext const&);
0136       void postGlobalEndLumi(GlobalContext const&);
0137 
0138       void preGlobalWriteLumi(GlobalContext const&);
0139       void postGlobalWriteLumi(GlobalContext const&);
0140 
0141       void preStreamBeginLumi(StreamContext const&);
0142       void postStreamBeginLumi(StreamContext const&);
0143 
0144       void preStreamEndLumi(StreamContext const&);
0145       void postStreamEndLumi(StreamContext const&);
0146 
0147       void preEvent(StreamContext const&);
0148       void postEvent(StreamContext const&);
0149 
0150       void prePathEvent(StreamContext const&, PathContext const&);
0151       void postPathEvent(StreamContext const&, PathContext const&, HLTPathStatus const&);
0152 
0153       void preModuleConstruction(ModuleDescription const& md);
0154       void postModuleConstruction(ModuleDescription const& md);
0155 
0156       void preModuleDestruction(ModuleDescription const& md);
0157       void postModuleDestruction(ModuleDescription const& md);
0158 
0159       void preModuleBeginJob(ModuleDescription const& md);
0160       void postModuleBeginJob(ModuleDescription const& md);
0161 
0162       void preModuleEndJob(ModuleDescription const& md);
0163       void postModuleEndJob(ModuleDescription const& md);
0164 
0165       void preModuleEventPrefetching(StreamContext const&, ModuleCallingContext const&);
0166       void postModuleEventPrefetching(StreamContext const&, ModuleCallingContext const&);
0167       void preModuleEvent(StreamContext const&, ModuleCallingContext const&);
0168       void postModuleEvent(StreamContext const&, ModuleCallingContext const&);
0169       void preModuleEventAcquire(StreamContext const&, ModuleCallingContext const&);
0170       void postModuleEventAcquire(StreamContext const&, ModuleCallingContext const&);
0171       void preModuleEventDelayedGet(StreamContext const&, ModuleCallingContext const&);
0172       void postModuleEventDelayedGet(StreamContext const&, ModuleCallingContext const&);
0173       void preEventReadFromSource(StreamContext const&, ModuleCallingContext const&);
0174       void postEventReadFromSource(StreamContext const&, ModuleCallingContext const&);
0175       void preModuleTransformPrefetching(StreamContext const&, ModuleCallingContext const&);
0176       void postModuleTransformPrefetching(StreamContext const&, ModuleCallingContext const&);
0177       void preModuleTransform(StreamContext const&, ModuleCallingContext const&);
0178       void postModuleTransform(StreamContext const&, ModuleCallingContext const&);
0179       void preModuleTransformAcquiring(StreamContext const&, ModuleCallingContext const&);
0180       void postModuleTransformAcquiring(StreamContext const&, ModuleCallingContext const&);
0181 
0182       void preModuleStreamPrefetching(StreamContext const&, ModuleCallingContext const&);
0183       void postModuleStreamPrefetching(StreamContext const&, ModuleCallingContext const&);
0184 
0185       void preModuleStreamBeginRun(StreamContext const&, ModuleCallingContext const&);
0186       void postModuleStreamBeginRun(StreamContext const&, ModuleCallingContext const&);
0187       void preModuleStreamEndRun(StreamContext const&, ModuleCallingContext const&);
0188       void postModuleStreamEndRun(StreamContext const&, ModuleCallingContext const&);
0189 
0190       void preModuleStreamBeginLumi(StreamContext const&, ModuleCallingContext const&);
0191       void postModuleStreamBeginLumi(StreamContext const&, ModuleCallingContext const&);
0192       void preModuleStreamEndLumi(StreamContext const&, ModuleCallingContext const&);
0193       void postModuleStreamEndLumi(StreamContext const&, ModuleCallingContext const&);
0194 
0195       void preModuleBeginProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0196       void postModuleBeginProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0197       void preModuleAccessInputProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0198       void postModuleAccessInputProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0199       void preModuleEndProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0200       void postModuleEndProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0201 
0202       void preModuleGlobalPrefetching(GlobalContext const&, ModuleCallingContext const&);
0203       void postModuleGlobalPrefetching(GlobalContext const&, ModuleCallingContext const&);
0204 
0205       void preModuleGlobalBeginRun(GlobalContext const&, ModuleCallingContext const&);
0206       void postModuleGlobalBeginRun(GlobalContext const&, ModuleCallingContext const&);
0207       void preModuleGlobalEndRun(GlobalContext const&, ModuleCallingContext const&);
0208       void postModuleGlobalEndRun(GlobalContext const&, ModuleCallingContext const&);
0209 
0210       void preModuleGlobalBeginLumi(GlobalContext const&, ModuleCallingContext const&);
0211       void postModuleGlobalBeginLumi(GlobalContext const&, ModuleCallingContext const&);
0212       void preModuleGlobalEndLumi(GlobalContext const&, ModuleCallingContext const&);
0213       void postModuleGlobalEndLumi(GlobalContext const&, ModuleCallingContext const&);
0214 
0215       void preModuleWriteProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0216       void postModuleWriteProcessBlock(GlobalContext const&, ModuleCallingContext const&);
0217 
0218       void preModuleWriteRun(GlobalContext const&, ModuleCallingContext const&);
0219       void postModuleWriteRun(GlobalContext const&, ModuleCallingContext const&);
0220 
0221       void preModuleWriteLumi(GlobalContext const&, ModuleCallingContext const&);
0222       void postModuleWriteLumi(GlobalContext const&, ModuleCallingContext const&);
0223 
0224       void preSourceConstruction(ModuleDescription const& md);
0225       void postSourceConstruction(ModuleDescription const& md);
0226 
0227       void preESModulePrefetching(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0228       void postESModulePrefetching(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0229       void preESModule(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0230       void postESModule(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0231       void preESModuleAcquire(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0232       void postESModuleAcquire(eventsetup::EventSetupRecordKey const&, ESModuleCallingContext const&);
0233 
0234       void printPaths(LogAbsolute& out, PathsAndConsumesOfModulesBase const& pathsAndConsumes) const;
0235 
0236       void printEDModulesWithDependentModules(LogAbsolute& out,
0237                                               std::vector<ModuleDescription const*> const& allModules,
0238                                               PathsAndConsumesOfModulesBase const& pathsAndConsumes) const;
0239 
0240       void printESModulesWithDependentModules(LogAbsolute& out,
0241                                               std::vector<eventsetup::ComponentDescription const*> const& allESModules,
0242                                               PathsAndConsumesOfModulesBase const& pathsAndConsumes) const;
0243 
0244       void printEDModulesWithConsumes(LogAbsolute& out,
0245                                       std::vector<ModuleDescription const*> const& allModules,
0246                                       PathsAndConsumesOfModulesBase const& pathsAndConsumes) const;
0247 
0248       void printESModulesWithConsumes(LogAbsolute& out,
0249                                       std::vector<eventsetup::ComponentDescription const*> const& allESModules,
0250                                       PathsAndConsumesOfModulesBase const& pathsAndConsumes) const;
0251 
0252     private:
0253       std::string indention_;
0254       std::set<std::string> dumpContextForLabels_;
0255       bool dumpNonModuleContext_;
0256       bool dumpPathsAndConsumes_;
0257       bool printTimestamps_;
0258       bool dumpEventSetupInfo_;
0259     };
0260   }  // namespace service
0261 }  // namespace edm
0262 
0263 using namespace edm::service;
0264 
0265 namespace {
0266 
0267   class TimeStamper {
0268   public:
0269     TimeStamper(bool enable) : enabled_(enable) {}
0270 
0271     friend std::ostream& operator<<(std::ostream& out, TimeStamper const& timestamp) {
0272       if (timestamp.enabled_)
0273         out << std::setprecision(2) << edm::TimeOfDay() << "  ";
0274       return out;
0275     }
0276 
0277   private:
0278     bool enabled_;
0279   };
0280 
0281 }  // namespace
0282 
0283 Tracer::Tracer(ParameterSet const& iPS, ActivityRegistry& iRegistry)
0284     : indention_(iPS.getUntrackedParameter<std::string>("indention")),
0285       dumpContextForLabels_(),
0286       dumpNonModuleContext_(iPS.getUntrackedParameter<bool>("dumpNonModuleContext")),
0287       dumpPathsAndConsumes_(iPS.getUntrackedParameter<bool>("dumpPathsAndConsumes")),
0288       printTimestamps_(iPS.getUntrackedParameter<bool>("printTimestamps")),
0289       dumpEventSetupInfo_(iPS.getUntrackedParameter<bool>("dumpEventSetupInfo")) {
0290   tracer::setupFile(iPS.getUntrackedParameter<std::string>("fileName"), iRegistry);
0291 
0292   if (not iPS.getUntrackedParameter<bool>("useMessageLogger"))
0293     return;
0294   for (std::string& label : iPS.getUntrackedParameter<std::vector<std::string>>("dumpContextForLabels"))
0295     dumpContextForLabels_.insert(std::move(label));
0296 
0297   iRegistry.watchPreallocate(this, &Tracer::preallocate);
0298 
0299   iRegistry.watchPreBeginJob(this, &Tracer::preBeginJob);
0300   iRegistry.watchPostBeginJob(this, &Tracer::postBeginJob);
0301   iRegistry.watchPreEndJob(this, &Tracer::preEndJob);
0302   iRegistry.watchPostEndJob(this, &Tracer::postEndJob);
0303   iRegistry.watchLookupInitializationComplete(this, &Tracer::lookupInitializationComplete);
0304 
0305   iRegistry.watchPreBeginStream(this, &Tracer::preBeginStream);
0306   iRegistry.watchPostBeginStream(this, &Tracer::postBeginStream);
0307   iRegistry.watchPreEndStream(this, &Tracer::preEndStream);
0308   iRegistry.watchPostEndStream(this, &Tracer::postEndStream);
0309 
0310   iRegistry.watchPreSourceEvent(this, &Tracer::preSourceEvent);
0311   iRegistry.watchPostSourceEvent(this, &Tracer::postSourceEvent);
0312 
0313   iRegistry.watchPreSourceLumi(this, &Tracer::preSourceLumi);
0314   iRegistry.watchPostSourceLumi(this, &Tracer::postSourceLumi);
0315 
0316   iRegistry.watchPreSourceRun(this, &Tracer::preSourceRun);
0317   iRegistry.watchPostSourceRun(this, &Tracer::postSourceRun);
0318 
0319   iRegistry.watchPreSourceProcessBlock(this, &Tracer::preSourceProcessBlock);
0320   iRegistry.watchPostSourceProcessBlock(this, &Tracer::postSourceProcessBlock);
0321 
0322   iRegistry.watchPreOpenFile(this, &Tracer::preOpenFile);
0323   iRegistry.watchPostOpenFile(this, &Tracer::postOpenFile);
0324 
0325   iRegistry.watchPreCloseFile(this, &Tracer::preCloseFile);
0326   iRegistry.watchPostCloseFile(this, &Tracer::postCloseFile);
0327 
0328   iRegistry.watchPreModuleBeginStream(this, &Tracer::preModuleBeginStream);
0329   iRegistry.watchPostModuleBeginStream(this, &Tracer::postModuleBeginStream);
0330 
0331   iRegistry.watchPreModuleEndStream(this, &Tracer::preModuleEndStream);
0332   iRegistry.watchPostModuleEndStream(this, &Tracer::postModuleEndStream);
0333 
0334   iRegistry.watchPreBeginProcessBlock(this, &Tracer::preBeginProcessBlock);
0335   iRegistry.watchPostBeginProcessBlock(this, &Tracer::postBeginProcessBlock);
0336 
0337   iRegistry.watchPreAccessInputProcessBlock(this, &Tracer::preAccessInputProcessBlock);
0338   iRegistry.watchPostAccessInputProcessBlock(this, &Tracer::postAccessInputProcessBlock);
0339 
0340   iRegistry.watchPreEndProcessBlock(this, &Tracer::preEndProcessBlock);
0341   iRegistry.watchPostEndProcessBlock(this, &Tracer::postEndProcessBlock);
0342 
0343   iRegistry.watchPreWriteProcessBlock(this, &Tracer::preWriteProcessBlock);
0344   iRegistry.watchPostWriteProcessBlock(this, &Tracer::postWriteProcessBlock);
0345 
0346   iRegistry.watchPreGlobalBeginRun(this, &Tracer::preGlobalBeginRun);
0347   iRegistry.watchPostGlobalBeginRun(this, &Tracer::postGlobalBeginRun);
0348 
0349   iRegistry.watchPreGlobalEndRun(this, &Tracer::preGlobalEndRun);
0350   iRegistry.watchPostGlobalEndRun(this, &Tracer::postGlobalEndRun);
0351 
0352   iRegistry.watchPreGlobalWriteRun(this, &Tracer::preGlobalWriteRun);
0353   iRegistry.watchPostGlobalWriteRun(this, &Tracer::postGlobalWriteRun);
0354 
0355   iRegistry.watchPreStreamBeginRun(this, &Tracer::preStreamBeginRun);
0356   iRegistry.watchPostStreamBeginRun(this, &Tracer::postStreamBeginRun);
0357 
0358   iRegistry.watchPreStreamEndRun(this, &Tracer::preStreamEndRun);
0359   iRegistry.watchPostStreamEndRun(this, &Tracer::postStreamEndRun);
0360 
0361   iRegistry.watchPreGlobalBeginLumi(this, &Tracer::preGlobalBeginLumi);
0362   iRegistry.watchPostGlobalBeginLumi(this, &Tracer::postGlobalBeginLumi);
0363 
0364   iRegistry.watchPreGlobalEndLumi(this, &Tracer::preGlobalEndLumi);
0365   iRegistry.watchPostGlobalEndLumi(this, &Tracer::postGlobalEndLumi);
0366 
0367   iRegistry.watchPreGlobalWriteLumi(this, &Tracer::preGlobalWriteLumi);
0368   iRegistry.watchPostGlobalWriteLumi(this, &Tracer::postGlobalWriteLumi);
0369 
0370   iRegistry.watchPreStreamBeginLumi(this, &Tracer::preStreamBeginLumi);
0371   iRegistry.watchPostStreamBeginLumi(this, &Tracer::postStreamBeginLumi);
0372 
0373   iRegistry.watchPreStreamEndLumi(this, &Tracer::preStreamEndLumi);
0374   iRegistry.watchPostStreamEndLumi(this, &Tracer::postStreamEndLumi);
0375 
0376   iRegistry.watchPreEvent(this, &Tracer::preEvent);
0377   iRegistry.watchPostEvent(this, &Tracer::postEvent);
0378 
0379   iRegistry.watchPrePathEvent(this, &Tracer::prePathEvent);
0380   iRegistry.watchPostPathEvent(this, &Tracer::postPathEvent);
0381 
0382   iRegistry.watchPreModuleConstruction(this, &Tracer::preModuleConstruction);
0383   iRegistry.watchPostModuleConstruction(this, &Tracer::postModuleConstruction);
0384 
0385   iRegistry.watchPreModuleDestruction(this, &Tracer::preModuleDestruction);
0386   iRegistry.watchPostModuleDestruction(this, &Tracer::postModuleDestruction);
0387 
0388   iRegistry.watchPreModuleBeginJob(this, &Tracer::preModuleBeginJob);
0389   iRegistry.watchPostModuleBeginJob(this, &Tracer::postModuleBeginJob);
0390 
0391   iRegistry.watchPreModuleEndJob(this, &Tracer::preModuleEndJob);
0392   iRegistry.watchPostModuleEndJob(this, &Tracer::postModuleEndJob);
0393 
0394   iRegistry.watchPreModuleEventPrefetching(this, &Tracer::preModuleEventPrefetching);
0395   iRegistry.watchPostModuleEventPrefetching(this, &Tracer::postModuleEventPrefetching);
0396   iRegistry.watchPreModuleEvent(this, &Tracer::preModuleEvent);
0397   iRegistry.watchPostModuleEvent(this, &Tracer::postModuleEvent);
0398   iRegistry.watchPreModuleEventAcquire(this, &Tracer::preModuleEventAcquire);
0399   iRegistry.watchPostModuleEventAcquire(this, &Tracer::postModuleEventAcquire);
0400   iRegistry.watchPreModuleEventDelayedGet(this, &Tracer::preModuleEventDelayedGet);
0401   iRegistry.watchPostModuleEventDelayedGet(this, &Tracer::postModuleEventDelayedGet);
0402   iRegistry.watchPreEventReadFromSource(this, &Tracer::preEventReadFromSource);
0403   iRegistry.watchPostEventReadFromSource(this, &Tracer::postEventReadFromSource);
0404   iRegistry.watchPreModuleTransformPrefetching(this, &Tracer::preModuleTransformPrefetching);
0405   iRegistry.watchPostModuleTransformPrefetching(this, &Tracer::postModuleTransformPrefetching);
0406   iRegistry.watchPreModuleTransform(this, &Tracer::preModuleTransform);
0407   iRegistry.watchPostModuleTransform(this, &Tracer::postModuleTransform);
0408   iRegistry.watchPreModuleTransformAcquiring(this, &Tracer::preModuleTransformAcquiring);
0409   iRegistry.watchPostModuleTransformAcquiring(this, &Tracer::postModuleTransformAcquiring);
0410 
0411   iRegistry.watchPreModuleStreamPrefetching(this, &Tracer::preModuleStreamPrefetching);
0412   iRegistry.watchPostModuleStreamPrefetching(this, &Tracer::postModuleStreamPrefetching);
0413 
0414   iRegistry.watchPreModuleStreamBeginRun(this, &Tracer::preModuleStreamBeginRun);
0415   iRegistry.watchPostModuleStreamBeginRun(this, &Tracer::postModuleStreamBeginRun);
0416   iRegistry.watchPreModuleStreamEndRun(this, &Tracer::preModuleStreamEndRun);
0417   iRegistry.watchPostModuleStreamEndRun(this, &Tracer::postModuleStreamEndRun);
0418 
0419   iRegistry.watchPreModuleStreamBeginLumi(this, &Tracer::preModuleStreamBeginLumi);
0420   iRegistry.watchPostModuleStreamBeginLumi(this, &Tracer::postModuleStreamBeginLumi);
0421   iRegistry.watchPreModuleStreamEndLumi(this, &Tracer::preModuleStreamEndLumi);
0422   iRegistry.watchPostModuleStreamEndLumi(this, &Tracer::postModuleStreamEndLumi);
0423 
0424   iRegistry.watchPreModuleBeginProcessBlock(this, &Tracer::preModuleBeginProcessBlock);
0425   iRegistry.watchPostModuleBeginProcessBlock(this, &Tracer::postModuleBeginProcessBlock);
0426   iRegistry.watchPreModuleAccessInputProcessBlock(this, &Tracer::preModuleAccessInputProcessBlock);
0427   iRegistry.watchPostModuleAccessInputProcessBlock(this, &Tracer::postModuleAccessInputProcessBlock);
0428   iRegistry.watchPreModuleEndProcessBlock(this, &Tracer::preModuleEndProcessBlock);
0429   iRegistry.watchPostModuleEndProcessBlock(this, &Tracer::postModuleEndProcessBlock);
0430 
0431   iRegistry.watchPreModuleGlobalPrefetching(this, &Tracer::preModuleGlobalPrefetching);
0432   iRegistry.watchPostModuleGlobalPrefetching(this, &Tracer::postModuleGlobalPrefetching);
0433 
0434   iRegistry.watchPreModuleGlobalBeginRun(this, &Tracer::preModuleGlobalBeginRun);
0435   iRegistry.watchPostModuleGlobalBeginRun(this, &Tracer::postModuleGlobalBeginRun);
0436   iRegistry.watchPreModuleGlobalEndRun(this, &Tracer::preModuleGlobalEndRun);
0437   iRegistry.watchPostModuleGlobalEndRun(this, &Tracer::postModuleGlobalEndRun);
0438 
0439   iRegistry.watchPreModuleGlobalBeginLumi(this, &Tracer::preModuleGlobalBeginLumi);
0440   iRegistry.watchPostModuleGlobalBeginLumi(this, &Tracer::postModuleGlobalBeginLumi);
0441   iRegistry.watchPreModuleGlobalEndLumi(this, &Tracer::preModuleGlobalEndLumi);
0442   iRegistry.watchPostModuleGlobalEndLumi(this, &Tracer::postModuleGlobalEndLumi);
0443 
0444   iRegistry.watchPreModuleWriteProcessBlock(this, &Tracer::preModuleWriteProcessBlock);
0445   iRegistry.watchPostModuleWriteProcessBlock(this, &Tracer::postModuleWriteProcessBlock);
0446 
0447   iRegistry.watchPreModuleWriteRun(this, &Tracer::preModuleWriteRun);
0448   iRegistry.watchPostModuleWriteRun(this, &Tracer::postModuleWriteRun);
0449 
0450   iRegistry.watchPreModuleWriteLumi(this, &Tracer::preModuleWriteLumi);
0451   iRegistry.watchPostModuleWriteLumi(this, &Tracer::postModuleWriteLumi);
0452 
0453   iRegistry.watchPreSourceConstruction(this, &Tracer::preSourceConstruction);
0454   iRegistry.watchPostSourceConstruction(this, &Tracer::postSourceConstruction);
0455 
0456   iRegistry.watchPreESModulePrefetching(this, &Tracer::preESModulePrefetching);
0457   iRegistry.watchPostESModulePrefetching(this, &Tracer::postESModulePrefetching);
0458   iRegistry.watchPreESModule(this, &Tracer::preESModule);
0459   iRegistry.watchPostESModule(this, &Tracer::postESModule);
0460   iRegistry.watchPreESModuleAcquire(this, &Tracer::preESModuleAcquire);
0461   iRegistry.watchPostESModuleAcquire(this, &Tracer::postESModuleAcquire);
0462 
0463   iRegistry.preSourceEarlyTerminationSignal_.connect([this](edm::TerminationOrigin iOrigin) {
0464     LogAbsolute out("Tracer");
0465     out << TimeStamper(printTimestamps_);
0466     out << indention_ << indention_ << " early termination before processing transition";
0467   });
0468   iRegistry.preStreamEarlyTerminationSignal_.connect(
0469       [this](edm::StreamContext const& iContext, edm::TerminationOrigin iOrigin) {
0470         LogAbsolute out("Tracer");
0471         out << TimeStamper(printTimestamps_);
0472         if (iContext.eventID().luminosityBlock() == 0) {
0473           out << indention_ << indention_ << " early termination of run: stream = " << iContext.streamID()
0474               << " run = " << iContext.eventID().run();
0475         } else {
0476           if (iContext.eventID().event() == 0) {
0477             out << indention_ << indention_ << " early termination of stream lumi: stream = " << iContext.streamID()
0478                 << " run = " << iContext.eventID().run() << " lumi = " << iContext.eventID().luminosityBlock();
0479           } else {
0480             out << indention_ << indention_ << " early termination of event: stream = " << iContext.streamID()
0481                 << " run = " << iContext.eventID().run() << " lumi = " << iContext.eventID().luminosityBlock()
0482                 << " event = " << iContext.eventID().event();
0483           }
0484         }
0485         out << " : time = " << iContext.timestamp().value();
0486       });
0487   iRegistry.preGlobalEarlyTerminationSignal_.connect(
0488       [this](edm::GlobalContext const& iContext, edm::TerminationOrigin iOrigin) {
0489         LogAbsolute out("Tracer");
0490         out << TimeStamper(printTimestamps_);
0491         if (iContext.luminosityBlockID().value() == 0) {
0492           out << indention_ << indention_ << " early termination of global run " << iContext.luminosityBlockID().run();
0493         } else {
0494           out << indention_ << indention_
0495               << " early termination of global lumi run = " << iContext.luminosityBlockID().run()
0496               << " lumi = " << iContext.luminosityBlockID().luminosityBlock();
0497         }
0498         out << " : time = " << iContext.timestamp().value();
0499       });
0500 
0501   iRegistry.esSyncIOVQueuingSignal_.connect([this](edm::IOVSyncValue const& iSync) {
0502     LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_
0503                           << " queuing: EventSetup synchronization " << iSync.eventID();
0504   });
0505   iRegistry.preESSyncIOVSignal_.connect([this](edm::IOVSyncValue const& iSync) {
0506     LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_
0507                           << " pre: EventSetup synchronizing " << iSync.eventID();
0508   });
0509   iRegistry.postESSyncIOVSignal_.connect([this](edm::IOVSyncValue const& iSync) {
0510     LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_
0511                           << " post: EventSetup synchronizing " << iSync.eventID();
0512   });
0513 }
0514 
0515 void Tracer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0516   edm::ParameterSetDescription desc;
0517   desc.addUntracked<std::string>("indention", "++")
0518       ->setComment("Prefix characters for output. The characters are repeated to form the indentation.");
0519   desc.addUntracked<std::vector<std::string>>("dumpContextForLabels", std::vector<std::string>{})
0520       ->setComment(
0521           "Prints context information to cout for the module transitions associated with these modules' labels");
0522   desc.addUntracked<bool>("dumpNonModuleContext", false)
0523       ->setComment("Prints context information to cout for the transitions not associated with any module label");
0524   desc.addUntracked<bool>("dumpPathsAndConsumes", false)
0525       ->setComment(
0526           "Prints information to cout about paths, endpaths, products consumed by modules and the dependencies between "
0527           "modules created by the products they consume");
0528   desc.addUntracked<bool>("printTimestamps", false)->setComment("Prints a time stamp for every transition");
0529   desc.addUntracked<bool>("dumpEventSetupInfo", false)
0530       ->setComment(
0531           "Prints info 3 times when an event setup cache is filled, before the lock, after the lock, and after "
0532           "filling");
0533   desc.addUntracked<bool>("useMessageLogger", true)
0534       ->setComment("If true, information is sent via the MessageLogger. Only use false if `fileName` is used.");
0535   desc.addUntracked<std::string>("fileName", "")
0536       ->setComment("If fileName is given, write a compact representation of tracer info to the file.");
0537   descriptions.add("Tracer", desc);
0538   descriptions.setComment(
0539       "This service prints each phase the framework is processing, e.g. constructing a module,running a module, etc.");
0540 }
0541 
0542 void Tracer::preallocate(service::SystemBounds const& bounds) {
0543   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_
0544                         << " preallocate: " << bounds.maxNumberOfConcurrentRuns() << " concurrent runs, "
0545                         << bounds.maxNumberOfConcurrentLuminosityBlocks() << " concurrent luminosity sections, "
0546                         << bounds.maxNumberOfStreams() << " streams";
0547 }
0548 
0549 void Tracer::preBeginJob(ProcessContext const& pc) {
0550   LogAbsolute out("Tracer");
0551   out << TimeStamper(printTimestamps_) << indention_ << " starting: begin job";
0552 }
0553 
0554 void Tracer::postBeginJob() {
0555   LogAbsolute out("Tracer");
0556   out << TimeStamper(printTimestamps_) << indention_ << " finished: begin job";
0557 }
0558 
0559 void Tracer::preEndJob() {
0560   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " starting: end job";
0561 }
0562 
0563 void Tracer::postEndJob() {
0564   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " finished: end job";
0565 }
0566 
0567 void Tracer::lookupInitializationComplete(PathsAndConsumesOfModulesBase const& pathsAndConsumes,
0568                                           ProcessContext const& pc) {
0569   if (dumpPathsAndConsumes_) {
0570     LogAbsolute out("Tracer");
0571     out << TimeStamper(printTimestamps_) << indention_ << " lookupInitializationComplete";
0572     out << "\n"
0573         << "Process name = " << pc.processName() << "\n";
0574     std::vector<ModuleDescription const*> const& allModules = pathsAndConsumes.allModules();
0575     std::vector<eventsetup::ComponentDescription const*> const& allESModules = pathsAndConsumes.allESModules();
0576 
0577     printPaths(out, pathsAndConsumes);
0578     printEDModulesWithDependentModules(out, allModules, pathsAndConsumes);
0579     printESModulesWithDependentModules(out, allESModules, pathsAndConsumes);
0580     printEDModulesWithConsumes(out, allModules, pathsAndConsumes);
0581     printESModulesWithConsumes(out, allESModules, pathsAndConsumes);
0582   }
0583 }
0584 
0585 void Tracer::printPaths(LogAbsolute& out, PathsAndConsumesOfModulesBase const& pathsAndConsumes) const {
0586   out << "paths:\n";
0587   std::vector<std::string> const& paths = pathsAndConsumes.paths();
0588   if (paths.empty()) {
0589     out << "  --- there are no paths ---\n";
0590   }
0591   for (auto const& path : paths) {
0592     out << "  " << path << "\n";
0593   }
0594   out << "end paths:\n";
0595   std::vector<std::string> const& endpaths = pathsAndConsumes.endPaths();
0596   if (endpaths.empty()) {
0597     out << "  --- there are no endpaths ---\n";
0598   }
0599   for (auto const& endpath : endpaths) {
0600     out << "  " << endpath << "\n";
0601   }
0602   for (unsigned int j = 0; j < paths.size(); ++j) {
0603     std::vector<ModuleDescription const*> const& modulesOnPath = pathsAndConsumes.modulesOnPath(j);
0604     out << "modules on path " << paths.at(j) << ":\n";
0605     if (modulesOnPath.empty()) {
0606       out << "  --- there are no modules on this path ---\n";
0607     }
0608     for (auto const& desc : modulesOnPath) {
0609       out << "  " << desc->moduleLabel() << "\n";
0610     }
0611   }
0612   for (unsigned int j = 0; j < endpaths.size(); ++j) {
0613     std::vector<ModuleDescription const*> const& modulesOnEndPath = pathsAndConsumes.modulesOnEndPath(j);
0614     out << "modules on endpath " << endpaths.at(j) << ":\n";
0615     if (modulesOnEndPath.empty()) {
0616       out << "  --- there are no modules on this endpath ---\n";
0617     }
0618     for (auto const& desc : modulesOnEndPath) {
0619       out << "  " << desc->moduleLabel() << "\n";
0620     }
0621   }
0622 }
0623 
0624 void Tracer::printEDModulesWithDependentModules(LogAbsolute& out,
0625                                                 std::vector<ModuleDescription const*> const& allModules,
0626                                                 PathsAndConsumesOfModulesBase const& pathsAndConsumes) const {
0627   out << "All modules and modules in the current process whose products they consume:\n"
0628          "(This does not include modules from previous processes or the source)\n"
0629          "(Exclusively considers consumed Event and EventSetup products, not Run, Lumi, or ProcessBlock products)\n";
0630   if (allModules.empty()) {
0631     out << "  --- there are no modules ---\n";
0632   }
0633   for (auto const& module : allModules) {
0634     out << "  " << module->moduleName() << "/'" << module->moduleLabel() << "'\n";
0635     unsigned int moduleID = module->id();
0636     if (pathsAndConsumes.moduleDescription(moduleID) != module) {
0637       throw cms::Exception("Tracer")
0638           << "Tracer::printEDModulesWithDependentModules, moduleDescription returns incorrect value";
0639     }
0640     std::vector<ModuleDescription const*> const& modulesWhoseProductsAreConsumedBy =
0641         pathsAndConsumes.modulesWhoseProductsAreConsumedBy(moduleID);
0642     if (!modulesWhoseProductsAreConsumedBy.empty()) {
0643       out << "    consumes products from these modules:\n";
0644       for (auto const& producingModule : modulesWhoseProductsAreConsumedBy) {
0645         out << "      " << producingModule->moduleName() << "/'" << producingModule->moduleLabel() << "'\n";
0646       }
0647     }
0648 
0649     for (unsigned int i = 0; i < kNumberOfEventSetupTransitions; ++i) {
0650       Transition transition = static_cast<Transition>(i);
0651       std::vector<eventsetup::ComponentDescription const*> const& esModulesWhoseProductsAreConsumedBy =
0652           pathsAndConsumes.esModulesWhoseProductsAreConsumedBy(moduleID, transition);
0653       if (!esModulesWhoseProductsAreConsumedBy.empty()) {
0654         const char* transitionString = "Event";
0655         if (transition == Transition::BeginLuminosityBlock) {
0656           transitionString = "BeginLuminosityBlock";
0657         } else if (transition == Transition::EndLuminosityBlock) {
0658           transitionString = "EndLuminosityBlock";
0659         } else if (transition == Transition::BeginRun) {
0660           transitionString = "BeginRun";
0661         } else if (transition == Transition::EndRun) {
0662           transitionString = "EndRun";
0663         }
0664         out << "    consumes products during " << transitionString << " from these EventSetup modules:\n";
0665         for (auto const& component : esModulesWhoseProductsAreConsumedBy) {
0666           std::string moduleBaseType("ESProducer");
0667           if (component->isSource_) {
0668             moduleBaseType = "ESSource";
0669           } else if (component->isLooper_) {
0670             moduleBaseType = "Looper producing EventSetup data";
0671           }
0672           out << "      " << component->type_ << "/'" << component->label_ << "' " << moduleBaseType << "\n";
0673         }
0674       }
0675     }
0676   }
0677 }
0678 
0679 void Tracer::printESModulesWithDependentModules(LogAbsolute& out,
0680                                                 std::vector<eventsetup::ComponentDescription const*> const& allESModules,
0681                                                 PathsAndConsumesOfModulesBase const& pathsAndConsumes) const {
0682   out << "All EventSetup modules:\n";
0683   std::vector<std::vector<eventsetup::ComponentDescription const*>> const& esModules =
0684       pathsAndConsumes.esModulesWhoseProductsAreConsumedByESModule();
0685   auto it = esModules.begin();
0686   if (allESModules.empty()) {
0687     out << "  --- there are no EventSetup modules ---\n";
0688   }
0689   for (auto const& componentDescription : allESModules) {
0690     out << "  " << componentDescription->type_ << "/'"
0691         << (componentDescription->label_.empty() ? componentDescription->type_ : componentDescription->label_) << "'\n";
0692     if (!it->empty()) {
0693       out << "    consumes products from these EventSetup modules:\n";
0694       for (auto const& consumedComponentDescription : *it) {
0695         out << "      " << consumedComponentDescription->type_ << "/'"
0696             << (consumedComponentDescription->label_.empty() ? consumedComponentDescription->type_
0697                                                              : consumedComponentDescription->label_)
0698             << "'\n";
0699       }
0700     }
0701     ++it;
0702   }
0703 }
0704 
0705 void Tracer::printEDModulesWithConsumes(LogAbsolute& out,
0706                                         std::vector<ModuleDescription const*> const& allModules,
0707                                         PathsAndConsumesOfModulesBase const& pathsAndConsumes) const {
0708   out << "All modules (listed by class and label) and all their consumed products.\n"
0709          "Consumed products are listed by type, label, instance, process.\n"
0710          "For products not in the event, 'processBlock', 'run' or 'lumi' is added to indicate the TTree they "
0711          "are from.\n"
0712          "For products that are declared with mayConsume, 'may consume' is added.\n"
0713          "For products consumed for Views, 'element type' is added\n"
0714          "For products only read from previous processes, 'skip current process' is added\n"
0715          "Consumed EventSetup products are listed by record type, data type, product label,\n"
0716          "transition of request, whether module is ESSource, ESProducer, or ESProducerLooper,\n"
0717          "and the module type and module label of the EventSetup module producing the product.\n"
0718          "Next is the produceMethodID (counts setWhatProduced calls in ESProducers, otherwise 0).\n"
0719          "If the requested module label was specified, then there is a remark in parentheses\n"
0720          "stating whether it matched the actual module label of the producer.\n";
0721   if (allModules.empty()) {
0722     out << "  --- there are no modules ---\n";
0723   }
0724   for (auto const* module : allModules) {
0725     out << "  " << module->moduleName() << "/'" << module->moduleLabel() << "'\n";
0726     std::vector<ModuleConsumesInfo> consumesInfos = pathsAndConsumes.moduleConsumesInfos(module->id());
0727     if (!consumesInfos.empty()) {
0728       out << "    consumes:\n";
0729       for (auto const& info : consumesInfos) {
0730         out << "      " << info.type() << " '" << info.label() << "' '" << info.instance();
0731         out << "' '" << info.process() << "'";
0732         if (info.branchType() == InLumi) {
0733           out << ", lumi";
0734         } else if (info.branchType() == InRun) {
0735           out << ", run";
0736         } else if (info.branchType() == InProcess) {
0737           out << ", processBlock";
0738         }
0739         if (!info.alwaysGets()) {
0740           out << ", may consume";
0741         }
0742         if (info.kindOfType() == ELEMENT_TYPE) {
0743           out << ", element type";
0744         }
0745         if (info.skipCurrentProcess()) {
0746           out << ", skip current process";
0747         }
0748         out << "\n";
0749       }
0750     }
0751     std::vector<ModuleConsumesESInfo> moduleConsumesESInfos = pathsAndConsumes.moduleConsumesESInfos(module->id());
0752     if (!moduleConsumesESInfos.empty()) {
0753       out << "    consumes from EventSetup:\n";
0754       for (auto const& info : moduleConsumesESInfos) {
0755         out << "      " << info.eventSetupRecordType_ << " " << info.productType_ << " '" << info.productLabel_ << "'";
0756         const char* transitionString = "Event";
0757         Transition transition = info.transitionOfConsumer_;
0758         if (transition == Transition::BeginLuminosityBlock) {
0759           transitionString = "BeginLuminosityBlock";
0760         } else if (transition == Transition::EndLuminosityBlock) {
0761           transitionString = "EndLuminosityBlock";
0762         } else if (transition == Transition::BeginRun) {
0763           transitionString = "BeginRun";
0764         } else if (transition == Transition::EndRun) {
0765           transitionString = "EndRun";
0766         }
0767         out << " " << transitionString;
0768         if (info.moduleLabelMismatch_) {
0769           out << ", EventSetup module label mismatch, requested '" << info.requestedModuleLabel_ << "'";
0770         } else if (info.moduleType_.empty()) {
0771           out << ", No EventSetup module configured to produce this data.";
0772         } else {
0773           out << " " << info.moduleBaseType() << " " << info.moduleType_ << "/'" << info.moduleLabel_ << "' "
0774               << info.produceMethodIDOfProducer_;
0775           if (!info.requestedModuleLabel_.empty()) {
0776             assert(info.requestedModuleLabel_ == info.moduleLabel_);
0777             out << " (module label matches requested label)";
0778           }
0779         }
0780         out << "\n";
0781       }
0782     }
0783   }
0784 }
0785 
0786 void Tracer::printESModulesWithConsumes(LogAbsolute& out,
0787                                         std::vector<eventsetup::ComponentDescription const*> const& allESModules,
0788                                         PathsAndConsumesOfModulesBase const& pathsAndConsumes) const {
0789   out << "All EventSetup modules (listed by class and label) and all their consumed products.\n"
0790          "These modules can only consume EventSetup products. These are listed by record type,\n"
0791          "data type, product label, whether the module is ESSource, ESProducer, or ESProducerLooper,\n"
0792          "and the module type and module label of the EventSetup module producing the product.\n"
0793          "Next is the produceMethodID (counts setWhatProduced calls in ESProducers, otherwise 0).\n"
0794          "If the requested module label was specified, then there is a remark in parentheses\n"
0795          "stating whether it matched the actual module label of the producer.\n";
0796   if (allESModules.empty()) {
0797     out << "  --- there are no EventSetup modules ---\n";
0798   }
0799   for (auto const& componentDescription : allESModules) {
0800     out << "  " << componentDescription->type_ << "/'"
0801         << (componentDescription->label_.empty() ? componentDescription->type_ : componentDescription->label_) << "'\n";
0802     std::vector<std::vector<ESModuleConsumesInfo>> esModuleConsumesInfos =
0803         pathsAndConsumes.esModuleConsumesInfos(componentDescription->id_);
0804 
0805     if (!esModuleConsumesInfos.empty()) {
0806       unsigned int produceMethod = 0;
0807       for (auto const& infosForOneProduceMethod : esModuleConsumesInfos) {
0808         out << "    Consumed products for consumer produce method: " << produceMethod << "\n";
0809         if (infosForOneProduceMethod.empty()) {
0810           out << "      --- no products consumed ---\n";
0811         }
0812         for (auto const& info : infosForOneProduceMethod) {
0813           if (info.mayConsumes_) {
0814             if (info.mayConsumesFirstEntry_) {
0815               out << "      " << info.eventSetupRecordType_ << " " << info.productType_ << "\n"
0816                   << "        May Consumes, available products that might be consumed (they match record type and "
0817                      "product type):\n";
0818             }
0819             if (info.mayConsumesNoProducts_) {
0820               out << "          --- no available products that match ---\n";
0821             } else {
0822               out << "          '" << info.productLabel_ << "' " << info.moduleBaseType() << " " << info.moduleType_
0823                   << "/'" << info.moduleLabel_ << "' " << info.produceMethodIDOfProducer_ << "\n";
0824             }
0825           } else {
0826             out << "      " << info.eventSetupRecordType_ << " " << info.productType_ << " '" << info.productLabel_
0827                 << "'";
0828             if (info.moduleLabelMismatch_) {
0829               out << ", EventSetup module label mismatch, requested '" << info.requestedModuleLabel_ << "'";
0830             } else if (info.moduleType_.empty()) {
0831               out << ", No EventSetup module configured to produce this data.";
0832             } else {
0833               out << " " << info.moduleBaseType() << " " << info.moduleType_ << "/'" << info.moduleLabel_ << "' "
0834                   << info.produceMethodIDOfProducer_;
0835               if (!info.requestedModuleLabel_.empty()) {
0836                 assert(info.requestedModuleLabel_ == info.moduleLabel_);
0837                 out << " (module label matches requested label)";
0838               }
0839             }
0840             out << "\n";
0841           }
0842         }
0843         ++produceMethod;
0844       }
0845     }
0846   }
0847 }
0848 
0849 void Tracer::preBeginStream(StreamContext const& sc) {
0850   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " starting: begin stream " << sc.streamID();
0851 }
0852 
0853 void Tracer::postBeginStream(StreamContext const& sc) {
0854   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " finished: begin stream " << sc.streamID();
0855 }
0856 
0857 void Tracer::preEndStream(StreamContext const& sc) {
0858   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " starting: end stream " << sc.streamID();
0859 }
0860 
0861 void Tracer::postEndStream(StreamContext const& sc) {
0862   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << " finished: end stream " << sc.streamID();
0863 }
0864 
0865 void Tracer::preSourceEvent(StreamID sid) {
0866   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " starting: source event";
0867 }
0868 
0869 void Tracer::postSourceEvent(StreamID sid) {
0870   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " finished: source event";
0871 }
0872 
0873 void Tracer::preSourceLumi(LuminosityBlockIndex index) {
0874   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " starting: source lumi";
0875 }
0876 
0877 void Tracer::postSourceLumi(LuminosityBlockIndex index) {
0878   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " finished: source lumi";
0879 }
0880 
0881 void Tracer::preSourceRun(RunIndex index) {
0882   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " starting: source run";
0883 }
0884 
0885 void Tracer::postSourceRun(RunIndex index) {
0886   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_ << " finished: source run";
0887 }
0888 
0889 void Tracer::preSourceProcessBlock() {
0890   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_
0891                         << " starting: source process block";
0892 }
0893 
0894 void Tracer::postSourceProcessBlock(std::string const& processName) {
0895   LogAbsolute("Tracer") << TimeStamper(printTimestamps_) << indention_ << indention_
0896                         << " finished: source process block " << processName;
0897 }
0898 
0899 void Tracer::preOpenFile(std::string const& lfn) {
0900   LogAbsolute out("Tracer");
0901   out << TimeStamper(printTimestamps_);
0902   out << indention_ << indention_ << " starting: open input file: lfn = " << lfn;
0903 }
0904 
0905 void Tracer::postOpenFile(std::string const& lfn) {
0906   LogAbsolute out("Tracer");
0907   out << TimeStamper(printTimestamps_);
0908   out << indention_ << indention_ << " finished: open input file: lfn = " << lfn;
0909 }
0910 
0911 void Tracer::preCloseFile(std::string const& lfn) {
0912   LogAbsolute out("Tracer");
0913   out << TimeStamper(printTimestamps_);
0914   out << indention_ << indention_ << " starting: close input file: lfn = " << lfn;
0915 }
0916 void Tracer::postCloseFile(std::string const& lfn) {
0917   LogAbsolute out("Tracer");
0918   out << TimeStamper(printTimestamps_);
0919   out << indention_ << indention_ << " finished: close input file: lfn = " << lfn;
0920 }
0921 
0922 void Tracer::preModuleBeginStream(StreamContext const& sc, ModuleCallingContext const& mcc) {
0923   ModuleDescription const& desc = *mcc.moduleDescription();
0924   LogAbsolute out("Tracer");
0925   out << TimeStamper(printTimestamps_);
0926   out << indention_ << indention_ << " starting: begin stream for module: stream = " << sc.streamID() << " label = '"
0927       << desc.moduleLabel() << "' id = " << desc.id();
0928   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
0929     out << "\n" << sc;
0930     out << mcc;
0931   }
0932 }
0933 
0934 void Tracer::postModuleBeginStream(StreamContext const& sc, ModuleCallingContext const& mcc) {
0935   ModuleDescription const& desc = *mcc.moduleDescription();
0936   LogAbsolute out("Tracer");
0937   out << TimeStamper(printTimestamps_);
0938   out << indention_ << indention_ << " finished: begin stream for module: stream = " << sc.streamID() << " label = '"
0939       << desc.moduleLabel() << "' id = " << desc.id();
0940   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
0941     out << "\n" << sc;
0942     out << mcc;
0943   }
0944 }
0945 
0946 void Tracer::preModuleEndStream(StreamContext const& sc, ModuleCallingContext const& mcc) {
0947   ModuleDescription const& desc = *mcc.moduleDescription();
0948   LogAbsolute out("Tracer");
0949   out << TimeStamper(printTimestamps_);
0950   out << indention_ << indention_ << " starting: end stream for module: stream = " << sc.streamID() << " label = '"
0951       << desc.moduleLabel() << "' id = " << desc.id();
0952   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
0953     out << "\n" << sc;
0954     out << mcc;
0955   }
0956 }
0957 
0958 void Tracer::postModuleEndStream(StreamContext const& sc, ModuleCallingContext const& mcc) {
0959   ModuleDescription const& desc = *mcc.moduleDescription();
0960   LogAbsolute out("Tracer");
0961   out << TimeStamper(printTimestamps_);
0962   out << indention_ << indention_ << " finished: end stream for module: stream = " << sc.streamID() << " label = '"
0963       << desc.moduleLabel() << "' id = " << desc.id();
0964   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
0965     out << "\n" << sc;
0966     out << mcc;
0967   }
0968 }
0969 
0970 void Tracer::preBeginProcessBlock(GlobalContext const& gc) {
0971   LogAbsolute out("Tracer");
0972   out << indention_ << indention_ << " starting: begin process block";
0973   if (dumpNonModuleContext_) {
0974     out << "\n" << gc;
0975   }
0976 }
0977 
0978 void Tracer::postBeginProcessBlock(GlobalContext const& gc) {
0979   LogAbsolute out("Tracer");
0980   out << indention_ << indention_ << " finished: begin process block";
0981   if (dumpNonModuleContext_) {
0982     out << "\n" << gc;
0983   }
0984 }
0985 
0986 void Tracer::preAccessInputProcessBlock(GlobalContext const& gc) {
0987   LogAbsolute out("Tracer");
0988   out << indention_ << indention_ << " starting: access input process block";
0989   if (dumpNonModuleContext_) {
0990     out << "\n" << gc;
0991   }
0992 }
0993 
0994 void Tracer::postAccessInputProcessBlock(GlobalContext const& gc) {
0995   LogAbsolute out("Tracer");
0996   out << indention_ << indention_ << " finished: access input process block";
0997   if (dumpNonModuleContext_) {
0998     out << "\n" << gc;
0999   }
1000 }
1001 
1002 void Tracer::preEndProcessBlock(GlobalContext const& gc) {
1003   LogAbsolute out("Tracer");
1004   out << indention_ << indention_ << " starting: end process block";
1005   if (dumpNonModuleContext_) {
1006     out << "\n" << gc;
1007   }
1008 }
1009 
1010 void Tracer::postEndProcessBlock(GlobalContext const& gc) {
1011   LogAbsolute out("Tracer");
1012   out << indention_ << indention_ << " finished: end process block";
1013   if (dumpNonModuleContext_) {
1014     out << "\n" << gc;
1015   }
1016 }
1017 
1018 void Tracer::preWriteProcessBlock(GlobalContext const& gc) {
1019   LogAbsolute out("Tracer");
1020   out << indention_ << indention_ << " starting: write process block";
1021   if (dumpNonModuleContext_) {
1022     out << "\n" << gc;
1023   }
1024 }
1025 
1026 void Tracer::postWriteProcessBlock(GlobalContext const& gc) {
1027   LogAbsolute out("Tracer");
1028   out << indention_ << indention_ << " finished: write process block";
1029   if (dumpNonModuleContext_) {
1030     out << "\n" << gc;
1031   }
1032 }
1033 
1034 void Tracer::preGlobalBeginRun(GlobalContext const& gc) {
1035   LogAbsolute out("Tracer");
1036   out << TimeStamper(printTimestamps_);
1037   out << indention_ << indention_ << " starting: global begin run " << gc.luminosityBlockID().run()
1038       << " : time = " << gc.timestamp().value();
1039   if (dumpNonModuleContext_) {
1040     out << "\n" << gc;
1041   }
1042 }
1043 
1044 void Tracer::postGlobalBeginRun(GlobalContext const& gc) {
1045   LogAbsolute out("Tracer");
1046   out << TimeStamper(printTimestamps_);
1047   out << indention_ << indention_ << " finished: global begin run " << gc.luminosityBlockID().run()
1048       << " : time = " << gc.timestamp().value();
1049   if (dumpNonModuleContext_) {
1050     out << "\n" << gc;
1051   }
1052 }
1053 
1054 void Tracer::preGlobalEndRun(GlobalContext const& gc) {
1055   LogAbsolute out("Tracer");
1056   out << TimeStamper(printTimestamps_);
1057   out << indention_ << indention_ << " starting: global end run " << gc.luminosityBlockID().run()
1058       << " : time = " << gc.timestamp().value();
1059   if (dumpNonModuleContext_) {
1060     out << "\n" << gc;
1061   }
1062 }
1063 
1064 void Tracer::postGlobalEndRun(GlobalContext const& gc) {
1065   LogAbsolute out("Tracer");
1066   out << TimeStamper(printTimestamps_);
1067   out << indention_ << indention_ << " finished: global end run " << gc.luminosityBlockID().run()
1068       << " : time = " << gc.timestamp().value();
1069   if (dumpNonModuleContext_) {
1070     out << "\n" << gc;
1071   }
1072 }
1073 
1074 void Tracer::preGlobalWriteRun(GlobalContext const& gc) {
1075   LogAbsolute out("Tracer");
1076   out << TimeStamper(printTimestamps_);
1077   out << indention_ << indention_ << " starting: global write run " << gc.luminosityBlockID().run()
1078       << " : time = " << gc.timestamp().value();
1079   if (dumpNonModuleContext_) {
1080     out << "\n" << gc;
1081   }
1082 }
1083 
1084 void Tracer::postGlobalWriteRun(GlobalContext const& gc) {
1085   LogAbsolute out("Tracer");
1086   out << TimeStamper(printTimestamps_);
1087   out << indention_ << indention_ << " finished: global write run " << gc.luminosityBlockID().run()
1088       << " : time = " << gc.timestamp().value();
1089   if (dumpNonModuleContext_) {
1090     out << "\n" << gc;
1091   }
1092 }
1093 
1094 void Tracer::preStreamBeginRun(StreamContext const& sc) {
1095   LogAbsolute out("Tracer");
1096   out << TimeStamper(printTimestamps_);
1097   out << indention_ << indention_ << " starting: begin run: stream = " << sc.streamID()
1098       << " run = " << sc.eventID().run() << " time = " << sc.timestamp().value();
1099   if (dumpNonModuleContext_) {
1100     out << "\n" << sc;
1101   }
1102 }
1103 
1104 void Tracer::postStreamBeginRun(StreamContext const& sc) {
1105   LogAbsolute out("Tracer");
1106   out << TimeStamper(printTimestamps_);
1107   out << indention_ << indention_ << " finished: begin run: stream = " << sc.streamID()
1108       << " run = " << sc.eventID().run() << " time = " << sc.timestamp().value();
1109   if (dumpNonModuleContext_) {
1110     out << "\n" << sc;
1111   }
1112 }
1113 
1114 void Tracer::preStreamEndRun(StreamContext const& sc) {
1115   LogAbsolute out("Tracer");
1116   out << TimeStamper(printTimestamps_);
1117   out << indention_ << indention_ << " starting: end run: stream = " << sc.streamID() << " run = " << sc.eventID().run()
1118       << " time = " << sc.timestamp().value();
1119   if (dumpNonModuleContext_) {
1120     out << "\n" << sc;
1121   }
1122 }
1123 
1124 void Tracer::postStreamEndRun(StreamContext const& sc) {
1125   LogAbsolute out("Tracer");
1126   out << TimeStamper(printTimestamps_);
1127   out << indention_ << indention_ << " finished: end run: stream = " << sc.streamID() << " run = " << sc.eventID().run()
1128       << " time = " << sc.timestamp().value();
1129   if (dumpNonModuleContext_) {
1130     out << "\n" << sc;
1131   }
1132 }
1133 
1134 void Tracer::preGlobalBeginLumi(GlobalContext const& gc) {
1135   LogAbsolute out("Tracer");
1136   out << TimeStamper(printTimestamps_);
1137   out << indention_ << indention_ << " starting: global begin lumi: run = " << gc.luminosityBlockID().run()
1138       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1139   if (dumpNonModuleContext_) {
1140     out << "\n" << gc;
1141   }
1142 }
1143 
1144 void Tracer::postGlobalBeginLumi(GlobalContext const& gc) {
1145   LogAbsolute out("Tracer");
1146   out << TimeStamper(printTimestamps_);
1147   out << indention_ << indention_ << " finished: global begin lumi: run = " << gc.luminosityBlockID().run()
1148       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1149   if (dumpNonModuleContext_) {
1150     out << "\n" << gc;
1151   }
1152 }
1153 
1154 void Tracer::preGlobalEndLumi(GlobalContext const& gc) {
1155   LogAbsolute out("Tracer");
1156   out << TimeStamper(printTimestamps_);
1157   out << indention_ << indention_ << " starting: global end lumi: run = " << gc.luminosityBlockID().run()
1158       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1159   if (dumpNonModuleContext_) {
1160     out << "\n" << gc;
1161   }
1162 }
1163 
1164 void Tracer::postGlobalEndLumi(GlobalContext const& gc) {
1165   LogAbsolute out("Tracer");
1166   out << TimeStamper(printTimestamps_);
1167   out << indention_ << indention_ << " finished: global end lumi: run = " << gc.luminosityBlockID().run()
1168       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1169   if (dumpNonModuleContext_) {
1170     out << "\n" << gc;
1171   }
1172 }
1173 
1174 void Tracer::preGlobalWriteLumi(GlobalContext const& gc) {
1175   LogAbsolute out("Tracer");
1176   out << TimeStamper(printTimestamps_);
1177   out << indention_ << indention_ << " starting: global write lumi: run = " << gc.luminosityBlockID().run()
1178       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1179   if (dumpNonModuleContext_) {
1180     out << "\n" << gc;
1181   }
1182 }
1183 
1184 void Tracer::postGlobalWriteLumi(GlobalContext const& gc) {
1185   LogAbsolute out("Tracer");
1186   out << TimeStamper(printTimestamps_);
1187   out << indention_ << indention_ << " finished: global write lumi: run = " << gc.luminosityBlockID().run()
1188       << " lumi = " << gc.luminosityBlockID().luminosityBlock() << " time = " << gc.timestamp().value();
1189   if (dumpNonModuleContext_) {
1190     out << "\n" << gc;
1191   }
1192 }
1193 
1194 void Tracer::preStreamBeginLumi(StreamContext const& sc) {
1195   LogAbsolute out("Tracer");
1196   out << TimeStamper(printTimestamps_);
1197   out << indention_ << indention_ << " starting: begin lumi: stream = " << sc.streamID()
1198       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1199       << " time = " << sc.timestamp().value();
1200   if (dumpNonModuleContext_) {
1201     out << "\n" << sc;
1202   }
1203 }
1204 
1205 void Tracer::postStreamBeginLumi(StreamContext const& sc) {
1206   LogAbsolute out("Tracer");
1207   out << TimeStamper(printTimestamps_);
1208   out << indention_ << indention_ << " finished: begin lumi: stream = " << sc.streamID()
1209       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1210       << " time = " << sc.timestamp().value();
1211   if (dumpNonModuleContext_) {
1212     out << "\n" << sc;
1213   }
1214 }
1215 
1216 void Tracer::preStreamEndLumi(StreamContext const& sc) {
1217   LogAbsolute out("Tracer");
1218   out << TimeStamper(printTimestamps_);
1219   out << indention_ << indention_ << " starting: end lumi: stream = " << sc.streamID()
1220       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1221       << " time = " << sc.timestamp().value();
1222   if (dumpNonModuleContext_) {
1223     out << "\n" << sc;
1224   }
1225 }
1226 
1227 void Tracer::postStreamEndLumi(StreamContext const& sc) {
1228   LogAbsolute out("Tracer");
1229   out << TimeStamper(printTimestamps_);
1230   out << indention_ << indention_ << " finished: end lumi: stream = " << sc.streamID()
1231       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1232       << " time = " << sc.timestamp().value();
1233   if (dumpNonModuleContext_) {
1234     out << "\n" << sc;
1235   }
1236 }
1237 
1238 void Tracer::preEvent(StreamContext const& sc) {
1239   LogAbsolute out("Tracer");
1240   out << TimeStamper(printTimestamps_);
1241   out << indention_ << indention_ << " starting: processing event : stream = " << sc.streamID()
1242       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1243       << " event = " << sc.eventID().event() << " time = " << sc.timestamp().value();
1244   if (dumpNonModuleContext_) {
1245     out << "\n" << sc;
1246   }
1247 }
1248 
1249 void Tracer::postEvent(StreamContext const& sc) {
1250   LogAbsolute out("Tracer");
1251   out << TimeStamper(printTimestamps_);
1252   out << indention_ << indention_ << " finished: processing event : stream = " << sc.streamID()
1253       << " run = " << sc.eventID().run() << " lumi = " << sc.eventID().luminosityBlock()
1254       << " event = " << sc.eventID().event() << " time = " << sc.timestamp().value();
1255   if (dumpNonModuleContext_) {
1256     out << "\n" << sc;
1257   }
1258 }
1259 
1260 void Tracer::prePathEvent(StreamContext const& sc, PathContext const& pc) {
1261   LogAbsolute out("Tracer");
1262   out << TimeStamper(printTimestamps_);
1263   out << indention_ << indention_ << indention_ << " starting: processing path '" << pc.pathName()
1264       << "' : stream = " << sc.streamID();
1265   if (dumpNonModuleContext_) {
1266     out << "\n" << sc;
1267     out << pc;
1268   }
1269 }
1270 
1271 void Tracer::postPathEvent(StreamContext const& sc, PathContext const& pc, HLTPathStatus const& hlts) {
1272   LogAbsolute out("Tracer");
1273   out << TimeStamper(printTimestamps_);
1274   out << indention_ << indention_ << indention_ << " finished: processing path '" << pc.pathName()
1275       << "' : stream = " << sc.streamID();
1276   if (dumpNonModuleContext_) {
1277     out << "\n" << sc;
1278     out << pc;
1279   }
1280 }
1281 
1282 void Tracer::preModuleConstruction(ModuleDescription const& desc) {
1283   LogAbsolute out("Tracer");
1284   out << TimeStamper(printTimestamps_);
1285   out << indention_ << indention_ << " starting: constructing module with label '" << desc.moduleLabel()
1286       << "' id = " << desc.id();
1287   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1288     out << "\n" << desc;
1289   }
1290 }
1291 
1292 void Tracer::postModuleConstruction(ModuleDescription const& desc) {
1293   LogAbsolute out("Tracer");
1294   out << TimeStamper(printTimestamps_);
1295   out << indention_ << indention_ << " finished: constructing module with label '" << desc.moduleLabel()
1296       << "' id = " << desc.id();
1297   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1298     out << "\n" << desc;
1299   }
1300 }
1301 
1302 void Tracer::preModuleDestruction(ModuleDescription const& desc) {
1303   LogAbsolute out("Tracer");
1304   out << TimeStamper(printTimestamps_);
1305   out << indention_ << indention_ << " starting: destructing module with label '" << desc.moduleLabel()
1306       << "' id = " << desc.id();
1307   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1308     out << "\n" << desc;
1309   }
1310 }
1311 
1312 void Tracer::postModuleDestruction(ModuleDescription const& desc) {
1313   LogAbsolute out("Tracer");
1314   out << TimeStamper(printTimestamps_);
1315   out << indention_ << indention_ << " finished: destructing module with label '" << desc.moduleLabel()
1316       << "' id = " << desc.id();
1317   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1318     out << "\n" << desc;
1319   }
1320 }
1321 
1322 void Tracer::preModuleBeginJob(ModuleDescription const& desc) {
1323   LogAbsolute out("Tracer");
1324   out << TimeStamper(printTimestamps_);
1325   out << indention_ << indention_;
1326   out << " starting: begin job for module with label '" << desc.moduleLabel() << "' id = " << desc.id();
1327   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1328     out << "\n" << desc;
1329   }
1330 }
1331 
1332 void Tracer::postModuleBeginJob(ModuleDescription const& desc) {
1333   LogAbsolute out("Tracer");
1334   out << TimeStamper(printTimestamps_);
1335   out << indention_ << indention_;
1336   out << " finished: begin job for module with label '" << desc.moduleLabel() << "' id = " << desc.id();
1337   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1338     out << "\n" << desc;
1339   }
1340 }
1341 
1342 void Tracer::preModuleEndJob(ModuleDescription const& desc) {
1343   LogAbsolute out("Tracer");
1344   out << TimeStamper(printTimestamps_);
1345   out << indention_ << indention_;
1346   out << " starting: end job for module with label '" << desc.moduleLabel() << "' id = " << desc.id();
1347   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1348     out << "\n" << desc;
1349   }
1350 }
1351 
1352 void Tracer::postModuleEndJob(ModuleDescription const& desc) {
1353   LogAbsolute out("Tracer");
1354   out << TimeStamper(printTimestamps_);
1355   out << indention_ << indention_;
1356   out << " finished: end job for module with label '" << desc.moduleLabel() << "' id = " << desc.id();
1357   if (dumpContextForLabels_.find(desc.moduleLabel()) != dumpContextForLabels_.end()) {
1358     out << "\n" << desc;
1359   }
1360 }
1361 
1362 void Tracer::preModuleEventPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1363   LogAbsolute out("Tracer");
1364   out << TimeStamper(printTimestamps_);
1365   unsigned int nIndents = mcc.depth() + 4;
1366   for (unsigned int i = 0; i < nIndents; ++i) {
1367     out << indention_;
1368   }
1369   out << " starting: prefetching before processing event for module: stream = " << sc.streamID() << " label = '"
1370       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1371   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1372     out << "\n" << sc;
1373     out << mcc;
1374   }
1375 }
1376 
1377 void Tracer::postModuleEventPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1378   LogAbsolute out("Tracer");
1379   out << TimeStamper(printTimestamps_);
1380   unsigned int nIndents = mcc.depth() + 4;
1381   for (unsigned int i = 0; i < nIndents; ++i) {
1382     out << indention_;
1383   }
1384   out << " finished: prefetching before processing event for module: stream = " << sc.streamID() << " label = '"
1385       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1386   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1387     out << "\n" << sc;
1388     out << mcc;
1389   }
1390 }
1391 
1392 void Tracer::preModuleEvent(StreamContext const& sc, ModuleCallingContext const& mcc) {
1393   LogAbsolute out("Tracer");
1394   out << TimeStamper(printTimestamps_);
1395   unsigned int nIndents = mcc.depth() + 4;
1396   for (unsigned int i = 0; i < nIndents; ++i) {
1397     out << indention_;
1398   }
1399   out << " starting: processing event for module: stream = " << sc.streamID() << " label = '"
1400       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1401   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1402     out << "\n" << sc;
1403     out << mcc;
1404   }
1405 }
1406 
1407 void Tracer::postModuleEvent(StreamContext const& sc, ModuleCallingContext const& mcc) {
1408   LogAbsolute out("Tracer");
1409   out << TimeStamper(printTimestamps_);
1410   unsigned int nIndents = mcc.depth() + 4;
1411   for (unsigned int i = 0; i < nIndents; ++i) {
1412     out << indention_;
1413   }
1414   out << " finished: processing event for module: stream = " << sc.streamID() << " label = '"
1415       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1416   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1417     out << "\n" << sc;
1418     out << mcc;
1419   }
1420 }
1421 
1422 void Tracer::preModuleEventAcquire(StreamContext const& sc, ModuleCallingContext const& mcc) {
1423   LogAbsolute out("Tracer");
1424   out << TimeStamper(printTimestamps_);
1425   unsigned int nIndents = mcc.depth() + 4;
1426   for (unsigned int i = 0; i < nIndents; ++i) {
1427     out << indention_;
1428   }
1429   out << " starting: processing event acquire for module: stream = " << sc.streamID() << " label = '"
1430       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1431 }
1432 
1433 void Tracer::postModuleEventAcquire(StreamContext const& sc, ModuleCallingContext const& mcc) {
1434   LogAbsolute out("Tracer");
1435   out << TimeStamper(printTimestamps_);
1436   unsigned int nIndents = mcc.depth() + 4;
1437   for (unsigned int i = 0; i < nIndents; ++i) {
1438     out << indention_;
1439   }
1440   out << " finished: processing event acquire for module: stream = " << sc.streamID() << " label = '"
1441       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1442 }
1443 
1444 void Tracer::preModuleEventDelayedGet(StreamContext const& sc, ModuleCallingContext const& mcc) {
1445   LogAbsolute out("Tracer");
1446   out << TimeStamper(printTimestamps_);
1447   unsigned int nIndents = mcc.depth() + 4;
1448   for (unsigned int i = 0; i < nIndents; ++i) {
1449     out << indention_;
1450   }
1451   out << " starting: delayed get while processing event for module: stream = " << sc.streamID() << " label = '"
1452       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1453   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1454     out << "\n" << sc;
1455     out << mcc;
1456   }
1457 }
1458 
1459 void Tracer::postModuleEventDelayedGet(StreamContext const& sc, ModuleCallingContext const& mcc) {
1460   LogAbsolute out("Tracer");
1461   out << TimeStamper(printTimestamps_);
1462   unsigned int nIndents = mcc.depth() + 4;
1463   for (unsigned int i = 0; i < nIndents; ++i) {
1464     out << indention_;
1465   }
1466   out << " finished: delayed get while processing event for module: stream = " << sc.streamID() << " label = '"
1467       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1468   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1469     out << "\n" << sc;
1470     out << mcc;
1471   }
1472 }
1473 
1474 void Tracer::preEventReadFromSource(StreamContext const& sc, ModuleCallingContext const& mcc) {
1475   LogAbsolute out("Tracer");
1476   out << TimeStamper(printTimestamps_);
1477   unsigned int nIndents = mcc.depth() + 5;
1478   for (unsigned int i = 0; i < nIndents; ++i) {
1479     out << indention_;
1480   }
1481   out << " starting: event delayed read from source: stream = " << sc.streamID() << " label = '"
1482       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1483 }
1484 
1485 void Tracer::postEventReadFromSource(StreamContext const& sc, ModuleCallingContext const& mcc) {
1486   LogAbsolute out("Tracer");
1487   out << TimeStamper(printTimestamps_);
1488   unsigned int nIndents = mcc.depth() + 5;
1489   for (unsigned int i = 0; i < nIndents; ++i) {
1490     out << indention_;
1491   }
1492   out << " finished: event delayed read from source: stream = " << sc.streamID() << " label = '"
1493       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1494 }
1495 
1496 void Tracer::preModuleTransformPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1497   LogAbsolute out("Tracer");
1498   out << TimeStamper(printTimestamps_);
1499   unsigned int nIndents = mcc.depth() + 4;
1500   for (unsigned int i = 0; i < nIndents; ++i) {
1501     out << indention_;
1502   }
1503   out << " starting: prefetching before transform in event for module: stream = " << sc.streamID() << " label = '"
1504       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1505       << " callID = " << mcc.callID();
1506   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1507     out << "\n" << sc;
1508     out << mcc;
1509   }
1510 }
1511 
1512 void Tracer::postModuleTransformPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1513   LogAbsolute out("Tracer");
1514   out << TimeStamper(printTimestamps_);
1515   unsigned int nIndents = mcc.depth() + 4;
1516   for (unsigned int i = 0; i < nIndents; ++i) {
1517     out << indention_;
1518   }
1519   out << " finished: prefetching before transform in event for module: stream = " << sc.streamID() << " label = '"
1520       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1521       << " callID = " << mcc.callID();
1522   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1523     out << "\n" << sc;
1524     out << mcc;
1525   }
1526 }
1527 
1528 void Tracer::preModuleTransform(StreamContext const& sc, ModuleCallingContext const& mcc) {
1529   LogAbsolute out("Tracer");
1530   out << TimeStamper(printTimestamps_);
1531   unsigned int nIndents = mcc.depth() + 4;
1532   for (unsigned int i = 0; i < nIndents; ++i) {
1533     out << indention_;
1534   }
1535   out << " starting: transform in event for module: stream = " << sc.streamID() << " label = '"
1536       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1537       << " callID = " << mcc.callID();
1538   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1539     out << "\n" << sc;
1540     out << mcc;
1541   }
1542 }
1543 
1544 void Tracer::postModuleTransform(StreamContext const& sc, ModuleCallingContext const& mcc) {
1545   LogAbsolute out("Tracer");
1546   out << TimeStamper(printTimestamps_);
1547   unsigned int nIndents = mcc.depth() + 4;
1548   for (unsigned int i = 0; i < nIndents; ++i) {
1549     out << indention_;
1550   }
1551   out << " finished: transform in event for module: stream = " << sc.streamID() << " label = '"
1552       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1553       << " callID = " << mcc.callID();
1554   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1555     out << "\n" << sc;
1556     out << mcc;
1557   }
1558 }
1559 
1560 void Tracer::preModuleTransformAcquiring(StreamContext const& sc, ModuleCallingContext const& mcc) {
1561   LogAbsolute out("Tracer");
1562   out << TimeStamper(printTimestamps_);
1563   unsigned int nIndents = mcc.depth() + 4;
1564   for (unsigned int i = 0; i < nIndents; ++i) {
1565     out << indention_;
1566   }
1567   out << " starting: acquiring before transform in event for module: stream = " << sc.streamID() << " label = '"
1568       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1569       << " callID = " << mcc.callID();
1570 }
1571 
1572 void Tracer::postModuleTransformAcquiring(StreamContext const& sc, ModuleCallingContext const& mcc) {
1573   LogAbsolute out("Tracer");
1574   out << TimeStamper(printTimestamps_);
1575   unsigned int nIndents = mcc.depth() + 4;
1576   for (unsigned int i = 0; i < nIndents; ++i) {
1577     out << indention_;
1578   }
1579   out << " finished: acquiring before transform in event acquire for module: stream = " << sc.streamID() << " label = '"
1580       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id()
1581       << " callID = " << mcc.callID();
1582 }
1583 void Tracer::preModuleStreamPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1584   LogAbsolute out("Tracer");
1585   out << TimeStamper(printTimestamps_);
1586   unsigned int nIndents = mcc.depth() + 4;
1587   for (unsigned int i = 0; i < nIndents; ++i) {
1588     out << indention_;
1589   }
1590   out << " starting: prefetching before processing " << transitionName(sc.transition())
1591       << " for module: stream = " << sc.streamID() << " label = '" << mcc.moduleDescription()->moduleLabel()
1592       << "' id = " << mcc.moduleDescription()->id();
1593   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1594     out << "\n" << sc;
1595     out << mcc;
1596   }
1597 }
1598 
1599 void Tracer::postModuleStreamPrefetching(StreamContext const& sc, ModuleCallingContext const& mcc) {
1600   LogAbsolute out("Tracer");
1601   out << TimeStamper(printTimestamps_);
1602   unsigned int nIndents = mcc.depth() + 4;
1603   for (unsigned int i = 0; i < nIndents; ++i) {
1604     out << indention_;
1605   }
1606   out << " finished: prefetching before processing " << transitionName(sc.transition())
1607       << " for module: stream = " << sc.streamID() << " label = '" << mcc.moduleDescription()->moduleLabel()
1608       << "' id = " << mcc.moduleDescription()->id();
1609   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1610     out << "\n" << sc;
1611     out << mcc;
1612   }
1613 }
1614 
1615 void Tracer::preModuleStreamBeginRun(StreamContext const& sc, ModuleCallingContext const& mcc) {
1616   LogAbsolute out("Tracer");
1617   out << TimeStamper(printTimestamps_);
1618   unsigned int nIndents = mcc.depth() + 3;
1619   for (unsigned int i = 0; i < nIndents; ++i) {
1620     out << indention_;
1621   }
1622   out << " starting: begin run for module: stream = " << sc.streamID() << " label = '"
1623       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1624   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1625     out << "\n" << sc;
1626     out << mcc;
1627   }
1628 }
1629 
1630 void Tracer::postModuleStreamBeginRun(StreamContext const& sc, ModuleCallingContext const& mcc) {
1631   LogAbsolute out("Tracer");
1632   out << TimeStamper(printTimestamps_);
1633   unsigned int nIndents = mcc.depth() + 3;
1634   for (unsigned int i = 0; i < nIndents; ++i) {
1635     out << indention_;
1636   }
1637   out << " finished: begin run for module: stream = " << sc.streamID() << " label = '"
1638       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1639   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1640     out << "\n" << sc;
1641     out << mcc;
1642   }
1643 }
1644 
1645 void Tracer::preModuleStreamEndRun(StreamContext const& sc, ModuleCallingContext const& mcc) {
1646   LogAbsolute out("Tracer");
1647   out << TimeStamper(printTimestamps_);
1648   unsigned int nIndents = mcc.depth() + 3;
1649   for (unsigned int i = 0; i < nIndents; ++i) {
1650     out << indention_;
1651   }
1652   out << " starting: end run for module: stream = " << sc.streamID() << " label = '"
1653       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1654   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1655     out << "\n" << sc;
1656     out << mcc;
1657   }
1658 }
1659 
1660 void Tracer::postModuleStreamEndRun(StreamContext const& sc, ModuleCallingContext const& mcc) {
1661   LogAbsolute out("Tracer");
1662   out << TimeStamper(printTimestamps_);
1663   unsigned int nIndents = mcc.depth() + 3;
1664   for (unsigned int i = 0; i < nIndents; ++i) {
1665     out << indention_;
1666   }
1667   out << " finished: end run for module: stream = " << sc.streamID() << " label = '"
1668       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1669   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1670     out << "\n" << sc;
1671     out << mcc;
1672   }
1673 }
1674 
1675 void Tracer::preModuleStreamBeginLumi(StreamContext const& sc, ModuleCallingContext const& mcc) {
1676   LogAbsolute out("Tracer");
1677   out << TimeStamper(printTimestamps_);
1678   unsigned int nIndents = mcc.depth() + 3;
1679   for (unsigned int i = 0; i < nIndents; ++i) {
1680     out << indention_;
1681   }
1682   out << " starting: begin lumi for module: stream = " << sc.streamID() << " label = '"
1683       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1684   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1685     out << "\n" << sc;
1686     out << mcc;
1687   }
1688 }
1689 
1690 void Tracer::postModuleStreamBeginLumi(StreamContext const& sc, ModuleCallingContext const& mcc) {
1691   LogAbsolute out("Tracer");
1692   out << TimeStamper(printTimestamps_);
1693   unsigned int nIndents = mcc.depth() + 3;
1694   for (unsigned int i = 0; i < nIndents; ++i) {
1695     out << indention_;
1696   }
1697   out << " finished: begin lumi for module: stream = " << sc.streamID() << " label = '"
1698       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1699   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1700     out << "\n" << sc;
1701     out << mcc;
1702   }
1703 }
1704 
1705 void Tracer::preModuleStreamEndLumi(StreamContext const& sc, ModuleCallingContext const& mcc) {
1706   LogAbsolute out("Tracer");
1707   out << TimeStamper(printTimestamps_);
1708   unsigned int nIndents = mcc.depth() + 3;
1709   for (unsigned int i = 0; i < nIndents; ++i) {
1710     out << indention_;
1711   }
1712   out << " starting: end lumi for module: stream = " << sc.streamID() << " label = '"
1713       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1714   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1715     out << "\n" << sc;
1716     out << mcc;
1717   }
1718 }
1719 
1720 void Tracer::postModuleStreamEndLumi(StreamContext const& sc, ModuleCallingContext const& mcc) {
1721   LogAbsolute out("Tracer");
1722   out << TimeStamper(printTimestamps_);
1723   unsigned int nIndents = mcc.depth() + 3;
1724   for (unsigned int i = 0; i < nIndents; ++i) {
1725     out << indention_;
1726   }
1727   out << " finished: end lumi for module: stream = " << sc.streamID() << " label = '"
1728       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1729   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1730     out << "\n" << sc;
1731     out << mcc;
1732   }
1733 }
1734 
1735 void Tracer::preModuleGlobalPrefetching(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1736   LogAbsolute out("Tracer");
1737   out << TimeStamper(printTimestamps_);
1738   unsigned int nIndents = mcc.depth() + 4;
1739   for (unsigned int i = 0; i < nIndents; ++i) {
1740     out << indention_;
1741   }
1742   out << " starting: prefetching before processing " << transitionName(gc.transition()) << " for module: label = '"
1743       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1744   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1745     out << "\n" << gc;
1746     out << mcc;
1747   }
1748 }
1749 
1750 void Tracer::postModuleGlobalPrefetching(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1751   LogAbsolute out("Tracer");
1752   out << TimeStamper(printTimestamps_);
1753   unsigned int nIndents = mcc.depth() + 4;
1754   for (unsigned int i = 0; i < nIndents; ++i) {
1755     out << indention_;
1756   }
1757   out << " finished: prefetching before processing " << transitionName(gc.transition()) << " for module: label = '"
1758       << mcc.moduleDescription()->moduleLabel() << "' id = " << mcc.moduleDescription()->id();
1759   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1760     out << "\n" << gc;
1761     out << mcc;
1762   }
1763 }
1764 
1765 void Tracer::preModuleBeginProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1766   LogAbsolute out("Tracer");
1767   unsigned int nIndents = mcc.depth() + 3;
1768   for (unsigned int i = 0; i < nIndents; ++i) {
1769     out << indention_;
1770   }
1771   out << " starting: begin process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1772       << "' id = " << mcc.moduleDescription()->id();
1773   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1774     out << "\n" << gc;
1775     out << mcc;
1776   }
1777 }
1778 
1779 void Tracer::postModuleBeginProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1780   LogAbsolute out("Tracer");
1781   unsigned int nIndents = mcc.depth() + 3;
1782   for (unsigned int i = 0; i < nIndents; ++i) {
1783     out << indention_;
1784   }
1785   out << " finished: begin process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1786       << "' id = " << mcc.moduleDescription()->id();
1787   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1788     out << "\n" << gc;
1789     out << mcc;
1790   }
1791 }
1792 
1793 void Tracer::preModuleAccessInputProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1794   LogAbsolute out("Tracer");
1795   unsigned int nIndents = mcc.depth() + 3;
1796   for (unsigned int i = 0; i < nIndents; ++i) {
1797     out << indention_;
1798   }
1799   out << " starting: access input process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1800       << "' id = " << mcc.moduleDescription()->id();
1801   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1802     out << "\n" << gc;
1803     out << mcc;
1804   }
1805 }
1806 
1807 void Tracer::postModuleAccessInputProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1808   LogAbsolute out("Tracer");
1809   unsigned int nIndents = mcc.depth() + 3;
1810   for (unsigned int i = 0; i < nIndents; ++i) {
1811     out << indention_;
1812   }
1813   out << " finished: access input process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1814       << "' id = " << mcc.moduleDescription()->id();
1815   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1816     out << "\n" << gc;
1817     out << mcc;
1818   }
1819 }
1820 
1821 void Tracer::preModuleEndProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1822   LogAbsolute out("Tracer");
1823   unsigned int nIndents = mcc.depth() + 3;
1824   for (unsigned int i = 0; i < nIndents; ++i) {
1825     out << indention_;
1826   }
1827   out << " starting: end process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1828       << "' id = " << mcc.moduleDescription()->id();
1829   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1830     out << "\n" << gc;
1831     out << mcc;
1832   }
1833 }
1834 
1835 void Tracer::postModuleEndProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1836   LogAbsolute out("Tracer");
1837   unsigned int nIndents = mcc.depth() + 3;
1838   for (unsigned int i = 0; i < nIndents; ++i) {
1839     out << indention_;
1840   }
1841   out << " finished: end process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1842       << "' id = " << mcc.moduleDescription()->id();
1843   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1844     out << "\n" << gc;
1845     out << mcc;
1846   }
1847 }
1848 
1849 void Tracer::preModuleGlobalBeginRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1850   LogAbsolute out("Tracer");
1851   out << TimeStamper(printTimestamps_);
1852   unsigned int nIndents = mcc.depth() + 3;
1853   for (unsigned int i = 0; i < nIndents; ++i) {
1854     out << indention_;
1855   }
1856   out << " starting: global begin run for module: label = '" << mcc.moduleDescription()->moduleLabel()
1857       << "' id = " << mcc.moduleDescription()->id();
1858   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1859     out << "\n" << gc;
1860     out << mcc;
1861   }
1862 }
1863 
1864 void Tracer::postModuleGlobalBeginRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1865   LogAbsolute out("Tracer");
1866   out << TimeStamper(printTimestamps_);
1867   unsigned int nIndents = mcc.depth() + 3;
1868   for (unsigned int i = 0; i < nIndents; ++i) {
1869     out << indention_;
1870   }
1871   out << " finished: global begin run for module: label = '" << mcc.moduleDescription()->moduleLabel()
1872       << "' id = " << mcc.moduleDescription()->id();
1873   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1874     out << "\n" << gc;
1875     out << mcc;
1876   }
1877 }
1878 
1879 void Tracer::preModuleGlobalEndRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1880   LogAbsolute out("Tracer");
1881   out << TimeStamper(printTimestamps_);
1882   unsigned int nIndents = mcc.depth() + 3;
1883   for (unsigned int i = 0; i < nIndents; ++i) {
1884     out << indention_;
1885   }
1886   out << " starting: global end run for module: label = '" << mcc.moduleDescription()->moduleLabel()
1887       << "' id = " << mcc.moduleDescription()->id();
1888   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1889     out << "\n" << gc;
1890     out << mcc;
1891   }
1892 }
1893 
1894 void Tracer::postModuleGlobalEndRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1895   LogAbsolute out("Tracer");
1896   out << TimeStamper(printTimestamps_);
1897   unsigned int nIndents = mcc.depth() + 3;
1898   for (unsigned int i = 0; i < nIndents; ++i) {
1899     out << indention_;
1900   }
1901   out << " finished: global end run for module: label = '" << mcc.moduleDescription()->moduleLabel()
1902       << "' id = " << mcc.moduleDescription()->id();
1903   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1904     out << "\n" << gc;
1905     out << mcc;
1906   }
1907 }
1908 
1909 void Tracer::preModuleGlobalBeginLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1910   LogAbsolute out("Tracer");
1911   out << TimeStamper(printTimestamps_);
1912   unsigned int nIndents = mcc.depth() + 3;
1913   for (unsigned int i = 0; i < nIndents; ++i) {
1914     out << indention_;
1915   }
1916   out << " starting: global begin lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
1917       << "' id = " << mcc.moduleDescription()->id();
1918   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1919     out << "\n" << gc;
1920     out << mcc;
1921   }
1922 }
1923 
1924 void Tracer::postModuleGlobalBeginLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1925   LogAbsolute out("Tracer");
1926   out << TimeStamper(printTimestamps_);
1927   unsigned int nIndents = mcc.depth() + 3;
1928   for (unsigned int i = 0; i < nIndents; ++i) {
1929     out << indention_;
1930   }
1931   out << " finished: global begin lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
1932       << "' id = " << mcc.moduleDescription()->id();
1933   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1934     out << "\n" << gc;
1935     out << mcc;
1936   }
1937 }
1938 
1939 void Tracer::preModuleGlobalEndLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1940   LogAbsolute out("Tracer");
1941   out << TimeStamper(printTimestamps_);
1942   unsigned int nIndents = mcc.depth() + 3;
1943   for (unsigned int i = 0; i < nIndents; ++i) {
1944     out << indention_;
1945   }
1946   out << " starting: global end lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
1947       << "' id = " << mcc.moduleDescription()->id();
1948   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1949     out << "\n" << gc;
1950     out << mcc;
1951   }
1952 }
1953 
1954 void Tracer::postModuleGlobalEndLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1955   LogAbsolute out("Tracer");
1956   out << TimeStamper(printTimestamps_);
1957   unsigned int nIndents = mcc.depth() + 3;
1958   for (unsigned int i = 0; i < nIndents; ++i) {
1959     out << indention_;
1960   }
1961   out << " finished: global end lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
1962       << "' id = " << mcc.moduleDescription()->id();
1963   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1964     out << "\n" << gc;
1965     out << mcc;
1966   }
1967 }
1968 
1969 void Tracer::preModuleWriteProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1970   LogAbsolute out("Tracer");
1971   unsigned int nIndents = mcc.depth() + 3;
1972   for (unsigned int i = 0; i < nIndents; ++i) {
1973     out << indention_;
1974   }
1975   out << " starting: write process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1976       << "' id = " << mcc.moduleDescription()->id();
1977   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1978     out << "\n" << gc;
1979     out << mcc;
1980   }
1981 }
1982 
1983 void Tracer::postModuleWriteProcessBlock(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1984   LogAbsolute out("Tracer");
1985   unsigned int nIndents = mcc.depth() + 3;
1986   for (unsigned int i = 0; i < nIndents; ++i) {
1987     out << indention_;
1988   }
1989   out << " finished: write process block for module: label = '" << mcc.moduleDescription()->moduleLabel()
1990       << "' id = " << mcc.moduleDescription()->id();
1991   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
1992     out << "\n" << gc;
1993     out << mcc;
1994   }
1995 }
1996 
1997 void Tracer::preModuleWriteRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
1998   LogAbsolute out("Tracer");
1999   out << TimeStamper(printTimestamps_);
2000   unsigned int nIndents = mcc.depth() + 3;
2001   for (unsigned int i = 0; i < nIndents; ++i) {
2002     out << indention_;
2003   }
2004   out << " starting: write run for module: label = '" << mcc.moduleDescription()->moduleLabel()
2005       << "' id = " << mcc.moduleDescription()->id();
2006   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
2007     out << "\n" << gc;
2008     out << mcc;
2009   }
2010 }
2011 
2012 void Tracer::postModuleWriteRun(GlobalContext const& gc, ModuleCallingContext const& mcc) {
2013   LogAbsolute out("Tracer");
2014   out << TimeStamper(printTimestamps_);
2015   unsigned int nIndents = mcc.depth() + 3;
2016   for (unsigned int i = 0; i < nIndents; ++i) {
2017     out << indention_;
2018   }
2019   out << " finished: write run for module: label = '" << mcc.moduleDescription()->moduleLabel()
2020       << "' id = " << mcc.moduleDescription()->id();
2021   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
2022     out << "\n" << gc;
2023     out << mcc;
2024   }
2025 }
2026 
2027 void Tracer::preModuleWriteLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
2028   LogAbsolute out("Tracer");
2029   out << TimeStamper(printTimestamps_);
2030   unsigned int nIndents = mcc.depth() + 3;
2031   for (unsigned int i = 0; i < nIndents; ++i) {
2032     out << indention_;
2033   }
2034   out << " starting: write lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
2035       << "' id = " << mcc.moduleDescription()->id();
2036   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
2037     out << "\n" << gc;
2038     out << mcc;
2039   }
2040 }
2041 
2042 void Tracer::postModuleWriteLumi(GlobalContext const& gc, ModuleCallingContext const& mcc) {
2043   LogAbsolute out("Tracer");
2044   out << TimeStamper(printTimestamps_);
2045   unsigned int nIndents = mcc.depth() + 3;
2046   for (unsigned int i = 0; i < nIndents; ++i) {
2047     out << indention_;
2048   }
2049   out << " finished: write lumi for module: label = '" << mcc.moduleDescription()->moduleLabel()
2050       << "' id = " << mcc.moduleDescription()->id();
2051   if (dumpContextForLabels_.find(mcc.moduleDescription()->moduleLabel()) != dumpContextForLabels_.end()) {
2052     out << "\n" << gc;
2053     out << mcc;
2054   }
2055 }
2056 
2057 void Tracer::preSourceConstruction(ModuleDescription const& desc) {
2058   LogAbsolute out("Tracer");
2059   out << TimeStamper(printTimestamps_);
2060   out << indention_;
2061   out << " starting: constructing source: " << desc.moduleName();
2062   if (dumpNonModuleContext_) {
2063     out << "\n" << desc;
2064   }
2065 }
2066 
2067 void Tracer::postSourceConstruction(ModuleDescription const& desc) {
2068   LogAbsolute out("Tracer");
2069   out << TimeStamper(printTimestamps_);
2070   out << indention_;
2071   out << " finished: constructing source: " << desc.moduleName();
2072   if (dumpNonModuleContext_) {
2073     out << "\n" << desc;
2074   }
2075 }
2076 
2077 void Tracer::preESModulePrefetching(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2078   LogAbsolute out("Tracer");
2079   out << TimeStamper(printTimestamps_);
2080   unsigned int nIndents = mcc.depth() + 4;
2081   for (unsigned int i = 0; i < nIndents; ++i) {
2082     out << indention_;
2083   }
2084   out << " starting: prefetching for esmodule: label = '" << mcc.componentDescription()->label_
2085       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2086 }
2087 
2088 void Tracer::postESModulePrefetching(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2089   LogAbsolute out("Tracer");
2090   out << TimeStamper(printTimestamps_);
2091   unsigned int nIndents = mcc.depth() + 4;
2092   for (unsigned int i = 0; i < nIndents; ++i) {
2093     out << indention_;
2094   }
2095   out << " finished: prefetching for esmodule: label = '" << mcc.componentDescription()->label_
2096       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2097 }
2098 
2099 void Tracer::preESModule(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2100   LogAbsolute out("Tracer");
2101   out << TimeStamper(printTimestamps_);
2102   unsigned int nIndents = mcc.depth() + 4;
2103   for (unsigned int i = 0; i < nIndents; ++i) {
2104     out << indention_;
2105   }
2106   out << " starting: processing esmodule: label = '" << mcc.componentDescription()->label_
2107       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2108 }
2109 
2110 void Tracer::postESModule(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2111   LogAbsolute out("Tracer");
2112   out << TimeStamper(printTimestamps_);
2113   unsigned int nIndents = mcc.depth() + 4;
2114   for (unsigned int i = 0; i < nIndents; ++i) {
2115     out << indention_;
2116   }
2117   out << " finished: processing esmodule: label = '" << mcc.componentDescription()->label_
2118       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2119 }
2120 
2121 void Tracer::preESModuleAcquire(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2122   LogAbsolute out("Tracer");
2123   out << TimeStamper(printTimestamps_);
2124   unsigned int nIndents = mcc.depth() + 4;
2125   for (unsigned int i = 0; i < nIndents; ++i) {
2126     out << indention_;
2127   }
2128   out << " starting: processing esmodule acquire: label = '" << mcc.componentDescription()->label_
2129       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2130 }
2131 
2132 void Tracer::postESModuleAcquire(eventsetup::EventSetupRecordKey const& iKey, ESModuleCallingContext const& mcc) {
2133   LogAbsolute out("Tracer");
2134   out << TimeStamper(printTimestamps_);
2135   unsigned int nIndents = mcc.depth() + 4;
2136   for (unsigned int i = 0; i < nIndents; ++i) {
2137     out << indention_;
2138   }
2139   out << " finished: processing esmodule acquire: label = '" << mcc.componentDescription()->label_
2140       << "' type = " << mcc.componentDescription()->type_ << " in record = " << iKey.name();
2141 }
2142 
2143 using edm::service::Tracer;
2144 DEFINE_FWK_SERVICE(Tracer);