File indexing completed on 2024-10-17 04:57:59
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0024 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0025 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
0026 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0027 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0028 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0029 #include "FWCore/Framework/interface/ESProducer.h"
0030 #include "FWCore/Framework/interface/EventSetupRecordIntervalFinder.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/Framework/interface/ModuleFactory.h"
0033 #include "FWCore/Framework/interface/SourceFactory.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0038 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0039 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0040
0041
0042
0043
0044
0045 class SiPixelFakeLorentzAngleESSource : public edm::ESProducer, public edm::EventSetupRecordIntervalFinder {
0046 public:
0047 SiPixelFakeLorentzAngleESSource(const edm::ParameterSet&);
0048 ~SiPixelFakeLorentzAngleESSource() override = default;
0049 virtual std::unique_ptr<SiPixelLorentzAngle> produce(const SiPixelLorentzAngleRcd&);
0050
0051 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0052
0053 protected:
0054 void setIntervalFor(const edm::eventsetup::EventSetupRecordKey&,
0055 const edm::IOVSyncValue&,
0056 edm::ValidityInterval&) override;
0057
0058 private:
0059 int HVgroup(int panel, int module);
0060
0061
0062 const edm::FileInPath fp_;
0063 const edm::FileInPath t_topo_fp_;
0064 const std::string myLabel_;
0065 typedef std::vector<edm::ParameterSet> Parameters;
0066 Parameters BPixParameters_;
0067 Parameters FPixParameters_;
0068 Parameters ModuleParameters_;
0069
0070 float bPixLorentzAnglePerTesla_;
0071 float fPixLorentzAnglePerTesla_;
0072 };
0073
0074
0075
0076
0077 SiPixelFakeLorentzAngleESSource::SiPixelFakeLorentzAngleESSource(const edm::ParameterSet& conf_)
0078 : fp_(conf_.getParameter<edm::FileInPath>("file")),
0079 t_topo_fp_(conf_.getParameter<edm::FileInPath>("topologyInput")),
0080 myLabel_(conf_.getParameter<std::string>("appendToDataLabel")),
0081 BPixParameters_(conf_.getParameter<Parameters>("BPixParameters")),
0082 FPixParameters_(conf_.getParameter<Parameters>("FPixParameters")),
0083 ModuleParameters_(conf_.getParameter<Parameters>("ModuleParameters")),
0084 bPixLorentzAnglePerTesla_((float)conf_.getUntrackedParameter<double>("bPixLorentzAnglePerTesla", -9999.)),
0085 fPixLorentzAnglePerTesla_((float)conf_.getUntrackedParameter<double>("fPixLorentzAnglePerTesla", -9999.)) {
0086 edm::LogInfo("SiPixelFakeLorentzAngleESSource::SiPixelFakeLorentzAngleESSource");
0087
0088 setWhatProduced(this);
0089 findingRecord<SiPixelLorentzAngleRcd>();
0090 }
0091
0092 std::unique_ptr<SiPixelLorentzAngle> SiPixelFakeLorentzAngleESSource::produce(const SiPixelLorentzAngleRcd&) {
0093 using namespace edm::es;
0094 unsigned int nmodules = 0;
0095 SiPixelLorentzAngle* obj = new SiPixelLorentzAngle();
0096 SiPixelDetInfoFileReader reader(fp_.fullPath());
0097 const std::vector<uint32_t>& DetIds = reader.getAllDetIds();
0098
0099 TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(t_topo_fp_.fullPath());
0100
0101
0102 for (const auto& detit : DetIds) {
0103 nmodules++;
0104 const DetId detid(detit);
0105 auto rawId = detid.rawId();
0106 int found = 0;
0107 int side = tTopo.side(detid);
0108
0109
0110 if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
0111 int layer = tTopo.pxbLayer(detid);
0112
0113 int ladder = tTopo.pxbLadder(detid);
0114
0115 int module = tTopo.pxbModule(detid);
0116 if (module < 5) {
0117 side = 1;
0118 } else {
0119 side = 2;
0120 }
0121
0122 LogDebug("SiPixelFakeLorentzAngleESSource") << " pixel barrel:"
0123 << " layer=" << layer << " ladder=" << ladder << " module=" << module
0124 << " rawId=" << rawId << " side " << side;
0125
0126
0127 if (bPixLorentzAnglePerTesla_ != -9999.) {
0128 edm::LogInfo("SiPixelFakeLorentzAngleESSource")
0129 << " LA = " << bPixLorentzAnglePerTesla_ << " common for all bpix" << std::endl;
0130 if (!obj->putLorentzAngle(detid.rawId(), bPixLorentzAnglePerTesla_))
0131 edm::LogError("SiPixelFakeLorentzAngleESSource")
0132 << "ERROR!: detid " << rawId << " already exists" << std::endl;
0133 } else {
0134
0135 for (const auto& it : ModuleParameters_) {
0136 if (it.getParameter<unsigned int>("rawid") == detid.rawId()) {
0137 float lorentzangle = (float)it.getParameter<double>("angle");
0138 if (!found) {
0139 obj->putLorentzAngle(detid.rawId(), lorentzangle);
0140 edm::LogInfo("SiPixelFakeLorentzAngleESSource")
0141 << " LA= " << lorentzangle << " individual value " << detid.rawId() << std::endl;
0142 found = 1;
0143 } else
0144 edm::LogError("SiPixelFakeLorentzAngleESSource") << "ERROR!: detid already exists" << std::endl;
0145 }
0146 }
0147
0148
0149 for (const auto& it : BPixParameters_) {
0150 if (it.exists("layer"))
0151 if (it.getParameter<int>("layer") != layer)
0152 continue;
0153 if (it.exists("ladder"))
0154 if (it.getParameter<int>("ladder") != ladder)
0155 continue;
0156 if (it.exists("module"))
0157 if (it.getParameter<int>("module") != module)
0158 continue;
0159 if (it.exists("side"))
0160 if (it.getParameter<int>("side") != side)
0161 continue;
0162 if (!found) {
0163 float lorentzangle = (float)it.getParameter<double>("angle");
0164 obj->putLorentzAngle(detid.rawId(), lorentzangle);
0165 edm::LogInfo("SiPixelFakeLorentzAngleESSource") << " LA= " << lorentzangle << std::endl;
0166 found = 2;
0167 } else if (found == 1) {
0168 edm::LogWarning("SiPixelFakeLorentzAngleESSource")
0169 << "detid already given in ModuleParameters, skipping ..." << std::endl;
0170 } else
0171 edm::LogError("SiPixelFakeLorentzAngleESSource")
0172 << " ERROR!: detid " << rawId << " already exists" << std::endl;
0173 }
0174 }
0175 } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
0176
0177
0178
0179
0180 PixelEndcapName pen(detid, &tTopo, true);
0181
0182
0183
0184
0185 int disk = pen.diskName();
0186 int blade = pen.bladeName();
0187 int panel = pen.pannelName();
0188 int ring = pen.ringName();
0189
0190 LogDebug("SiPixelFakeLorentzAngleESSource") << " pixel endcap:"
0191 << " side=" << side << " disk=" << disk << " blade =" << blade
0192 << " pannel=" << panel << " ring=" << ring << " rawId=" << rawId;
0193
0194
0195 if (fPixLorentzAnglePerTesla_ != -9999.) {
0196 edm::LogInfo("SiPixelFakeLorentzAngleESSource")
0197 << " LA = " << fPixLorentzAnglePerTesla_ << " common for all fpix" << std::endl;
0198 if (!obj->putLorentzAngle(detid.rawId(), fPixLorentzAnglePerTesla_))
0199 edm::LogError("SiPixelFakeLorentzAngleESSource") << " ERROR! detid already exists" << std::endl;
0200 } else {
0201
0202 for (const auto& it : ModuleParameters_) {
0203 if (it.getParameter<unsigned int>("rawid") == detid.rawId()) {
0204 float lorentzangle = (float)it.getParameter<double>("angle");
0205 if (!found) {
0206 obj->putLorentzAngle(detid.rawId(), lorentzangle);
0207 edm::LogInfo("SiPixelFakeLorentzAngleESSource")
0208 << " LA= " << lorentzangle << " individual value " << detid.rawId() << std::endl;
0209 found = 1;
0210 } else
0211 edm::LogError("SiPixelFakeLorentzAngleESSource")
0212 << "ERROR!: detid " << rawId << " already exists" << std::endl;
0213 }
0214 }
0215
0216
0217 for (const auto& it : FPixParameters_) {
0218 if (it.exists("side"))
0219 if (it.getParameter<int>("side") != side)
0220 continue;
0221 if (it.exists("disk"))
0222 if (it.getParameter<int>("disk") != disk)
0223 continue;
0224 if (it.exists("ring"))
0225 if (it.getParameter<int>("ring") != ring)
0226 continue;
0227 if (it.exists("blade"))
0228 if (it.getParameter<int>("blade") != blade)
0229 continue;
0230 if (it.exists("panel"))
0231 if (it.getParameter<int>("panel") != panel)
0232 continue;
0233 if (it.exists("HVgroup"))
0234 if (it.getParameter<int>("HVgroup") != HVgroup(panel, ring))
0235 continue;
0236 if (!found) {
0237 float lorentzangle = (float)it.getParameter<double>("angle");
0238 obj->putLorentzAngle(detid.rawId(), lorentzangle);
0239 edm::LogInfo("SiPixelFakeLorentzAngleESSource") << " LA= " << lorentzangle << std::endl;
0240 found = 2;
0241 } else if (found == 1) {
0242 edm::LogWarning("SiPixelFakeLorentzAngleESSource")
0243 << " detid " << rawId << " already given in ModuleParameters, skipping ..." << std::endl;
0244 } else
0245 edm::LogError("SiPixelFakeLorentzAngleESSource")
0246 << "ERROR!: detid" << rawId << "already exists" << std::endl;
0247 }
0248 }
0249 }
0250 }
0251
0252 edm::LogInfo("SiPixelFakeLorentzAngleESSource") << "Modules = " << nmodules << std::endl;
0253
0254 return std::unique_ptr<SiPixelLorentzAngle>(obj);
0255 }
0256
0257 void SiPixelFakeLorentzAngleESSource::setIntervalFor(const edm::eventsetup::EventSetupRecordKey&,
0258 const edm::IOVSyncValue& iosv,
0259 edm::ValidityInterval& oValidity) {
0260 edm::ValidityInterval infinity(iosv.beginOfTime(), iosv.endOfTime());
0261 oValidity = infinity;
0262 }
0263
0264 int SiPixelFakeLorentzAngleESSource::HVgroup(int panel, int module) {
0265 if (1 == panel && (1 == module || 2 == module)) {
0266 return 1;
0267 } else if (1 == panel && (3 == module || 4 == module)) {
0268 return 2;
0269 } else if (2 == panel && 1 == module) {
0270 return 1;
0271 } else if (2 == panel && (2 == module || 3 == module)) {
0272 return 2;
0273 } else {
0274 edm::LogError("SiPixelFakeLorentzAngleESSource")
0275 << " ERROR! in SiPixelFakeLorentzAngleESSource::HVgroup(...), panel = " << panel << ", module = " << module
0276 << std::endl;
0277 return 0;
0278 }
0279 }
0280
0281 void SiPixelFakeLorentzAngleESSource::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0282 edm::ParameterSetDescription desc;
0283 desc.setComment("ESSource to supply per-module SiPixelLorentzAngle payloads in the EventSetup");
0284
0285 desc.add<edm::FileInPath>(
0286 "file", edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt"))
0287 ->setComment("Tracker skimmed geometry");
0288 desc.add<edm::FileInPath>("topologyInput",
0289 edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml"))
0290 ->setComment("Tracker Topology");
0291
0292 desc.add<std::string>("appendToDataLabel", "")->setComment("label to which write the data");
0293 desc.addUntracked<double>("bPixLorentzAnglePerTesla", -9999.)->setComment("LA value for all BPix");
0294 desc.addUntracked<double>("fPixLorentzAnglePerTesla", -9999.)->setComment("LA value for all FPix");
0295
0296 edm::ParameterSetDescription desc_BPixParameters;
0297 desc_BPixParameters.addOptional<int>("layer");
0298 desc_BPixParameters.addOptional<int>("ladder");
0299 desc_BPixParameters.addOptional<int>("module");
0300 desc_BPixParameters.addOptional<int>("side");
0301 desc_BPixParameters.add<double>("angle");
0302 std::vector<edm::ParameterSet> default_BPixParametersToAdd;
0303 desc.addVPSet("BPixParameters", desc_BPixParameters, default_BPixParametersToAdd)
0304 ->setComment("LA values for given BPix regions");
0305
0306 edm::ParameterSetDescription desc_FPixParameters;
0307 desc_FPixParameters.addOptional<int>("side");
0308 desc_FPixParameters.addOptional<int>("disk");
0309 desc_FPixParameters.addOptional<int>("ring");
0310 desc_FPixParameters.addOptional<int>("blade");
0311 desc_FPixParameters.addOptional<int>("panel");
0312 desc_FPixParameters.addOptional<int>("HVgroup");
0313 desc_FPixParameters.add<double>("angle");
0314 std::vector<edm::ParameterSet> default_FPixParametersToAdd;
0315 desc.addVPSet("FPixParameters", desc_FPixParameters, default_FPixParametersToAdd)
0316 ->setComment("LA values for given FPix regions");
0317
0318 edm::ParameterSetDescription desc_ModuleParameters;
0319 desc_ModuleParameters.add<unsigned int>("rawid");
0320 desc_ModuleParameters.add<double>("angle");
0321 std::vector<edm::ParameterSet> default_ModuleParametersToAdd;
0322 desc.addVPSet("ModuleParameters", desc_ModuleParameters, default_ModuleParametersToAdd)
0323 ->setComment("LA values for given modules");
0324
0325 descriptions.addWithDefaultLabel(desc);
0326 }
0327
0328 DEFINE_FWK_EVENTSETUP_SOURCE(SiPixelFakeLorentzAngleESSource);