CMS 3D CMS Logo

SiPixelVCalDB.cc
Go to the documentation of this file.
1 // system includes
2 #include <fstream>
3 #include <iostream>
4 #include <map>
5 #include <memory>
6 #include <string>
7 
8 // user includes
26 
28 public:
29  explicit SiPixelVCalDB(const edm::ParameterSet& conf);
30  explicit SiPixelVCalDB();
31  ~SiPixelVCalDB() override;
32  void analyze(const edm::Event&, const edm::EventSetup&) override;
33 
34 private:
38  typedef std::vector<edm::ParameterSet> Parameters;
41 };
42 
43 using namespace std;
44 using namespace edm;
45 
47  : tkGeomToken_(esConsumes()), tkTopoToken_(esConsumes()) {
48  recordName_ = iConfig.getUntrackedParameter<std::string>("record", "SiPixelVCalRcd");
49  BPixParameters_ = iConfig.getUntrackedParameter<Parameters>("BPixParameters");
50  FPixParameters_ = iConfig.getUntrackedParameter<Parameters>("FPixParameters");
51 }
52 
54 
55 // Analyzer: Functions that gets called by framework every event
56 void SiPixelVCalDB::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
57  SiPixelVCal vcal;
58  bool phase1 = true;
59 
60  // Retrieve tracker topology from geometry
61  const TrackerTopology* const tTopo = &iSetup.getData(tkTopoToken_);
62 
63  // Retrieve old style tracker geometry from geometry
64  const TrackerGeometry* pDD = &iSetup.getData(tkGeomToken_);
65  edm::LogPrint("SiPixelVCalDB") << " There are " << pDD->detUnits().size() << " modules" << std::endl;
66 
67  for (const auto& it : pDD->detUnits()) {
68  if (dynamic_cast<PixelGeomDetUnit const*>(it) != nullptr) {
69  const DetId detid = it->geographicalId();
70  const unsigned int rawDetId = detid.rawId();
71  int subid = detid.subdetId();
72 
73  // FILL BPIX
74  if (subid == static_cast<int>(PixelSubdetector::PixelBarrel)) {
75  int layer = tTopo->pxbLayer(detid); // 1, 2, 3, 4
76  int ladder = tTopo->pxbLadder(detid); // 1-12/28/44/64
77  edm::LogPrint("SiPixelVCalDB") << " pixel barrel:"
78  << " detId=" << rawDetId << ", layer=" << layer << ", ladder=" << ladder;
79  for (Parameters::iterator it = BPixParameters_.begin(); it != BPixParameters_.end(); ++it) {
80  if (it->getParameter<int>("layer") == layer && it->getParameter<int>("ladder") == ladder) {
81  float slope = (float)it->getParameter<double>("slope");
82  float offset = (float)it->getParameter<double>("offset");
83  edm::LogPrint("SiPixelVCalDB") << "; VCal slope " << slope << ", offset " << offset;
84  // edm::LogInfo("SiPixelVCalDB") << " detId " << rawDetId << " \t
85  // VCal slope " << slope << ", offset " << offset;
86  vcal.putSlopeAndOffset(detid, slope, offset);
87  }
88  }
89  edm::LogPrint("SiPixelVCalDB") << std::endl;
90 
91  // FILL FPIX
92  } else if (subid == static_cast<int>(PixelSubdetector::PixelEndcap)) {
93  PixelEndcapName fpix(detid, tTopo, phase1);
94  int side = tTopo->pxfSide(detid); // 1 (-z), 2 for (+z)
95  int disk = fpix.diskName(); // 1, 2, 3
96  int disk2 = tTopo->pxfDisk(detid); // 1, 2, 3
97  int ring = fpix.ringName(); // 1 (lower), 2 (upper)
98  if (disk != disk2) {
99  edm::LogError("SiPixelVCalDB::analyze")
100  << "Found contradicting FPIX disk number: " << disk << " vs." << disk2 << std::endl;
101  }
102  edm::LogPrint("SiPixelVCalDB") << " pixel endcap:"
103  << " detId=" << rawDetId << ", side=" << side << ", disk=" << disk
104  << ", ring=" << ring;
105  for (Parameters::iterator it = FPixParameters_.begin(); it != FPixParameters_.end(); ++it) {
106  if (it->getParameter<int>("side") == side && it->getParameter<int>("disk") == disk &&
107  it->getParameter<int>("ring") == ring) {
108  float slope = (float)it->getParameter<double>("slope");
109  float offset = (float)it->getParameter<double>("offset");
110  edm::LogPrint("SiPixelVCalDB") << "; VCal slope " << slope << ", offset " << offset;
111  // edm::LogInfo("SiPixelVCalDB") << " detId " << rawDetId << " \t
112  // VCal slope " << slope << ", offset " << offset;
113  vcal.putSlopeAndOffset(rawDetId, slope, offset);
114  }
115  }
116  edm::LogPrint("SiPixelVCalDB") << std::endl;
117 
118  } else {
119  edm::LogError("SiPixelVCalDB::analyze") << "detid is Pixel but neither bpix nor fpix" << std::endl;
120  }
121  }
122  }
123 
124  // Save to DB
126  if (mydbservice.isAvailable()) {
127  try {
128  if (mydbservice->isNewTagRequest(recordName_)) {
129  mydbservice->createOneIOV<SiPixelVCal>(vcal, mydbservice->beginOfTime(), recordName_);
130  } else {
131  mydbservice->appendOneIOV<SiPixelVCal>(vcal, mydbservice->currentTime(), recordName_);
132  }
133  } catch (const cond::Exception& er) {
134  edm::LogError("SiPixelVCalDB") << er.what() << std::endl;
135  } catch (const std::exception& er) {
136  edm::LogError("SiPixelVCalDB") << "caught std::exception " << er.what() << std::endl;
137  } catch (...) {
138  edm::LogError("SiPixelVCalDB") << "Funny error" << std::endl;
139  }
140  } else {
141  edm::LogError("SiPixelVCalDB") << "Service is unavailable" << std::endl;
142  }
143 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int pxbLayer(const DetId &id) const
Base exception class for the object to relational access.
Definition: Exception.h:11
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
int ringName() const
ring Id
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
static const double slope[3]
Parameters FPixParameters_
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
unsigned int pxbLadder(const DetId &id) const
Log< level::Error, false > LogError
void createOneIOV(const T &payload, cond::Time_t firstSinceTime, const std::string &recordName)
T getUntrackedParameter(std::string const &, T const &) const
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
void appendOneIOV(const T &payload, cond::Time_t sinceTime, const std::string &recordName)
bool isNewTagRequest(const std::string &recordName)
int diskName() const
disk id
unsigned int pxfDisk(const DetId &id) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Log< level::Warning, true > LogPrint
Definition: DetId.h:17
unsigned int pxfSide(const DetId &id) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Parameters BPixParameters_
std::vector< edm::ParameterSet > Parameters
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
std::string recordName_
bool isAvailable() const
Definition: Service.h:40
char const * what() const noexcept override
Definition: Exception.cc:107
~SiPixelVCalDB() override
void putSlopeAndOffset(std::map< unsigned int, VCal > &vcal)
Definition: SiPixelVCal.h:22