File indexing completed on 2024-04-06 12:14:15
0001 #include <string>
0002 #include <vector>
0003
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Utilities/interface/ESGetToken.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "Utilities/General/interface/ClassName.h"
0010 #include "CondFormats/Alignment/interface/Alignments.h"
0011 #include "CondFormats/Alignment/interface/AlignTransform.h"
0012 #include "CondFormats/AlignmentRecord/interface/EBAlignmentRcd.h"
0013 #include "CondFormats/AlignmentRecord/interface/EEAlignmentRcd.h"
0014 #include "CondFormats/AlignmentRecord/interface/ESAlignmentRcd.h"
0015
0016 class CaloAlignmentRcdRead : public edm::one::EDAnalyzer<> {
0017 public:
0018 explicit CaloAlignmentRcdRead(const edm::ParameterSet& )
0019 : ebToken_{esConsumes<Alignments, EBAlignmentRcd>(edm::ESInputTag{})},
0020 eeToken_{esConsumes<Alignments, EEAlignmentRcd>(edm::ESInputTag{})},
0021 esToken_{esConsumes<Alignments, ESAlignmentRcd>(edm::ESInputTag{})},
0022 nEventCalls_(0) {}
0023 ~CaloAlignmentRcdRead() override {}
0024
0025 template <typename T>
0026 void dumpAlignments(const edm::EventSetup& evtSetup, edm::ESGetToken<Alignments, T>& token);
0027
0028 void beginJob() override {}
0029 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0030 void endJob() override {}
0031
0032 private:
0033 edm::ESGetToken<Alignments, EBAlignmentRcd> ebToken_;
0034 edm::ESGetToken<Alignments, EEAlignmentRcd> eeToken_;
0035 edm::ESGetToken<Alignments, ESAlignmentRcd> esToken_;
0036
0037 unsigned int nEventCalls_;
0038 };
0039
0040 template <typename T>
0041 void CaloAlignmentRcdRead::dumpAlignments(const edm::EventSetup& evtSetup, edm::ESGetToken<Alignments, T>& token) {
0042 const auto& alignments = evtSetup.getData(token);
0043
0044 std::string recordName = Demangle(typeid(T).name())();
0045
0046 LogDebug("CaloAlignmentRcdRead") << "Dumping alignments: " << recordName;
0047
0048 for (const auto& i : alignments.m_align) {
0049 LogDebug("CaloAlignmentRcdRead") << "entry " << i.rawId() << " translation " << i.translation() << " angles "
0050 << i.rotation().eulerAngles();
0051 }
0052 }
0053
0054 void CaloAlignmentRcdRead::analyze(const edm::Event& , const edm::EventSetup& evtSetup) {
0055 if (nEventCalls_ > 0) {
0056 edm::LogWarning("CaloAlignmentRcdRead") << "Reading from DB to be done only once, "
0057 << "set 'untracked PSet maxEvents = {untracked int32 input = 1}'.";
0058
0059 return;
0060 }
0061
0062 LogDebug("CaloAlignmentRcdRead") << "Reading from database in CaloAlignmentRcdRead::analyze...";
0063
0064 dumpAlignments(evtSetup, ebToken_);
0065 dumpAlignments(evtSetup, eeToken_);
0066 dumpAlignments(evtSetup, esToken_);
0067
0068 LogDebug("CaloAlignmentRcdRead") << "done!";
0069
0070 nEventCalls_++;
0071 }
0072
0073 DEFINE_FWK_MODULE(CaloAlignmentRcdRead);