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
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
|
/*
* See header file for a description of this class.
*
* Author of original version: M. Giunta
* \author A. Vilela Pereira
*/
#include "DTVDriftWriter.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "Geometry/DTGeometry/interface/DTGeometry.h"
#include "Geometry/Records/interface/MuonGeometryRecord.h"
#include "Geometry/DTGeometry/interface/DTSuperLayer.h"
#include "CondFormats/DTObjects/interface/DTMtime.h"
#include "CondFormats/DataRecord/interface/DTMtimeRcd.h"
#include "CondFormats/DTObjects/interface/DTRecoConditions.h"
#include "CondFormats/DataRecord/interface/DTRecoConditionsVdriftRcd.h"
#include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
#include "CalibMuon/DTCalibration/interface/DTVDriftPluginFactory.h"
#include "CalibMuon/DTCalibration/interface/DTVDriftBaseAlgo.h"
#include <string>
#include <vector>
using namespace std;
using namespace edm;
DTVDriftWriter::DTVDriftWriter(const ParameterSet& pset)
: mTimeMapToken_(esConsumes<edm::Transition::BeginRun>()),
vDriftMapToken_(esConsumes<edm::Transition::BeginRun>()),
dtGeomToken_(esConsumes<edm::Transition::BeginRun>()),
granularity_(pset.getUntrackedParameter<string>("calibGranularity", "bySL")),
mTimeMap_(nullptr),
vDriftMap_(nullptr),
vDriftAlgo_{DTVDriftPluginFactory::get()->create(pset.getParameter<string>("vDriftAlgo"),
pset.getParameter<ParameterSet>("vDriftAlgoConfig"),
consumesCollector())} {
LogVerbatim("Calibration") << "[DTVDriftWriter]Constructor called!";
if (granularity_ != "bySL")
throw cms::Exception("Configuration")
<< "[DTVDriftWriter] Check parameter calibGranularity: " << granularity_ << " option not available.";
readLegacyVDriftDB = pset.getParameter<bool>("readLegacyVDriftDB");
writeLegacyVDriftDB = pset.getParameter<bool>("writeLegacyVDriftDB");
}
DTVDriftWriter::~DTVDriftWriter() { LogVerbatim("Calibration") << "[DTVDriftWriter]Destructor called!"; }
void DTVDriftWriter::beginRun(const edm::Run& run, const edm::EventSetup& setup) {
// Get the map of vdrift from the Setup
if (readLegacyVDriftDB) {
mTimeMap_ = &setup.getData(mTimeMapToken_);
} else {
vDriftMap_ = &setup.getData(vDriftMapToken_);
// Consistency check: no parametrization is implemented for the time being
int version = vDriftMap_->version();
if (version != 1) {
throw cms::Exception("Configuration") << "only version 1 is presently supported for VDriftDB";
}
}
// Get geometry from Event Setup
dtGeom_ = setup.getHandle(dtGeomToken_);
// Pass EventSetup to concrete implementation
vDriftAlgo_->setES(setup);
}
void DTVDriftWriter::endJob() {
// Create the object to be written to DB
std::unique_ptr<DTMtime> mTimeNewMap;
std::unique_ptr<DTRecoConditions> vDriftNewMap;
if (writeLegacyVDriftDB) {
mTimeNewMap = std::make_unique<DTMtime>();
} else {
vDriftNewMap = std::make_unique<DTRecoConditions>();
vDriftNewMap->setFormulaExpr("[0]");
//vDriftNewMap->setFormulaExpr("[0]*(1-[1]*x)"); // add parametrization for dependency along Y
vDriftNewMap->setVersion(1);
}
if (granularity_ == "bySL") {
// Get all the sls from the geometry
const vector<const DTSuperLayer*>& superLayers = dtGeom_->superLayers();
auto sl = superLayers.begin();
auto sl_end = superLayers.end();
for (; sl != sl_end; ++sl) {
DTSuperLayerId slId = (*sl)->id();
// Compute vDrift
float vDriftNew = -1.;
float resolutionNew = -1;
try {
dtCalibration::DTVDriftData vDriftData = vDriftAlgo_->compute(slId);
vDriftNew = vDriftData.vdrift;
resolutionNew = vDriftData.resolution;
LogVerbatim("Calibration") << "vDrift for: " << slId << " Mean " << vDriftNew << " Resolution "
<< resolutionNew;
} catch (cms::Exception& e) { // Failure to compute new value, fall back to old table
LogError("Calibration") << e.explainSelf();
if (readLegacyVDriftDB) { //...reading old db format...
int status = mTimeMap_->get(slId, vDriftNew, resolutionNew, DTVelocityUnits::cm_per_ns);
if (status == 0) { // not found; silently skip this SL
continue;
}
} else { //...reading new db format
try {
vDriftNew = vDriftMap_->get(DTWireId(slId.rawId()));
} catch (cms::Exception& e2) {
// not found; silently skip this SL
continue;
}
}
LogVerbatim("Calibration") << "Keep original vDrift for: " << slId << " Mean " << vDriftNew << " Resolution "
<< resolutionNew;
}
// Add value to the vdrift table
if (writeLegacyVDriftDB) {
mTimeNewMap->set(slId, vDriftNew, resolutionNew, DTVelocityUnits::cm_per_ns);
} else {
vector<double> params = {vDriftNew};
vDriftNewMap->set(DTWireId(slId.rawId()), params);
}
} // End of loop on superlayers
}
// Write the vDrift object to DB
LogVerbatim("Calibration") << "[DTVDriftWriter]Writing vdrift object to DB!";
if (writeLegacyVDriftDB) {
string record = "DTMtimeRcd";
DTCalibDBUtils::writeToDB<DTMtime>(record, *mTimeNewMap);
} else {
DTCalibDBUtils::writeToDB<DTRecoConditions>("DTRecoConditionsVdriftRcd", *vDriftNewMap);
}
}
|