1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
|
/****************************************************************************
*
* Authors:
* Jan Kašpar (jan.kaspar@gmail.com)
*
****************************************************************************/
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Utilities/interface/Exception.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ESWatcher.h"
#include "CondFormats/AlignmentRecord/interface/RPRealAlignmentRecord.h"
#include "CondFormats/AlignmentRecord/interface/RPMisalignedAlignmentRecord.h"
#include "CondFormats/PPSObjects/interface/CTPPSRPAlignmentCorrectionsData.h"
//----------------------------------------------------------------------------------------------------
/**
* \brief Class to print out information on current geometry.
**/
class CTPPSAlignmentInfo : public edm::one::EDAnalyzer<> {
public:
explicit CTPPSAlignmentInfo(const edm::ParameterSet&);
private:
std::string alignmentType_;
edm::ESWatcher<RPRealAlignmentRecord> watcherRealAlignments_;
edm::ESWatcher<RPMisalignedAlignmentRecord> watcherMisalignedAlignments_;
edm::ESGetToken<CTPPSRPAlignmentCorrectionsData, RPRealAlignmentRecord> alignmentToken_;
void analyze(const edm::Event&, const edm::EventSetup&) override;
void printInfo(const CTPPSRPAlignmentCorrectionsData& alignments, const edm::Event& event) const;
};
CTPPSAlignmentInfo::CTPPSAlignmentInfo(const edm::ParameterSet& iConfig)
: alignmentType_(iConfig.getUntrackedParameter<std::string>("alignmentType", "real")),
alignmentToken_(esConsumes()) {}
//----------------------------------------------------------------------------------------------------
void CTPPSAlignmentInfo::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
if (alignmentType_ == "real") {
if (watcherRealAlignments_.check(iSetup)) {
auto const& alignments = iSetup.getData(alignmentToken_);
printInfo(alignments, iEvent);
}
return;
}
else if (alignmentType_ == "misaligned") {
if (watcherMisalignedAlignments_.check(iSetup)) {
auto const& alignments = iSetup.getData(alignmentToken_);
printInfo(alignments, iEvent);
}
return;
}
throw cms::Exception("CTPPSAlignmentInfo") << "Unknown geometry type: `" << alignmentType_ << "'.";
}
//----------------------------------------------------------------------------------------------------
void CTPPSAlignmentInfo::printInfo(const CTPPSRPAlignmentCorrectionsData& alignments, const edm::Event& event) const {
time_t unixTime = event.time().unixTime();
char timeStr[50];
strftime(timeStr, 50, "%F %T", localtime(&unixTime));
edm::LogInfo("CTPPSAlignmentInfo") << "New " << alignmentType_ << " alignments found in run=" << event.id().run()
<< ", event=" << event.id().event() << ", UNIX timestamp=" << unixTime << " ("
<< timeStr << "):\n"
<< alignments;
}
//----------------------------------------------------------------------------------------------------
DEFINE_FWK_MODULE(CTPPSAlignmentInfo);
|