CMS 3D CMS Logo

SiPixelLorentzAngleDB.cc
Go to the documentation of this file.
1 // system includes
2 #include <map>
3 #include <memory>
4 #include <string>
5 #include <iostream>
6 #include <fstream>
7 
8 // user includes
27 
29 public:
30  explicit SiPixelLorentzAngleDB(const edm::ParameterSet& conf);
31  ~SiPixelLorentzAngleDB() override;
32  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
33 
34 private:
37 
38  unsigned int HVgroup(unsigned int panel, unsigned int module);
39 
40  std::vector<std::pair<uint32_t, float> > detid_la;
43 
44  typedef std::vector<edm::ParameterSet> Parameters;
48 
50  bool useFile_;
51 };
52 
53 using namespace std;
54 using namespace edm;
55 
56 //Constructor
57 
59  : tkGeomToken_(esConsumes()), tkTopoToken_(esConsumes()) {
60  magneticField_ = conf.getParameter<double>("magneticField");
61  recordName_ = conf.getUntrackedParameter<std::string>("record", "SiPixelLorentzAngleRcd");
62  useFile_ = conf.getParameter<bool>("useFile");
63  fileName_ = conf.getParameter<string>("fileName");
64 
65  BPixParameters_ = conf.getUntrackedParameter<Parameters>("BPixParameters");
66  FPixParameters_ = conf.getUntrackedParameter<Parameters>("FPixParameters");
67  ModuleParameters_ = conf.getUntrackedParameter<Parameters>("ModuleParameters");
68 }
69 
70 // Virtual destructor needed.
72 
73 // Analyzer: Functions that gets called by framework every event
74 
77 
78  //Retrieve tracker topology from geometry
79  const TrackerTopology* tTopo = &es.getData(tkTopoToken_);
80 
81  //Retrieve old style tracker geometry from geometry
82  const TrackerGeometry* pDD = &es.getData(tkGeomToken_);
83  edm::LogInfo("SiPixelLorentzAngle (old)")
84  << " There are " << pDD->detUnits().size() << " detectors (old)" << std::endl;
85 
86  for (const auto& it : pDD->detUnits()) {
87  if (dynamic_cast<PixelGeomDetUnit const*>(it) != nullptr) {
88  DetId detid = it->geographicalId();
89  const DetId detidc = it->geographicalId();
90 
91  // fill bpix values for LA
92  if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
93  edm::LogPrint("SiPixelLorentzAngleDB")
94  << " pixel barrel:"
95  << " layer=" << tTopo->pxbLayer(detidc.rawId()) << " ladder=" << tTopo->pxbLadder(detidc.rawId())
96  << " module=" << tTopo->pxbModule(detidc.rawId()) << " rawId=" << detidc.rawId() << endl;
97 
98  if (!useFile_) {
99  //first individuals are put
100  for (Parameters::iterator it = ModuleParameters_.begin(); it != ModuleParameters_.end(); ++it) {
101  if (it->getParameter<unsigned int>("rawid") == detidc.rawId()) {
102  float lorentzangle = (float)it->getParameter<double>("angle");
103  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
104  edm::LogPrint("SiPixelLorentzAngleDB")
105  << " individual value=" << lorentzangle << " put into rawid=" << detid.rawId() << endl;
106  }
107  }
108 
109  //modules already put are automatically skipped
110  for (Parameters::iterator it = BPixParameters_.begin(); it != BPixParameters_.end(); ++it) {
111  if (it->getParameter<unsigned int>("module") == tTopo->pxbModule(detidc.rawId()) &&
112  it->getParameter<unsigned int>("layer") == tTopo->pxbLayer(detidc.rawId())) {
113  float lorentzangle = (float)it->getParameter<double>("angle");
114  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
115  }
116  }
117 
118  } else {
119  edm::LogError("SiPixelLorentzAngleDB")
120  << "[SiPixelLorentzAngleDB::analyze] method for reading file not implemented yet" << std::endl;
121  }
122 
123  // fill fpix values for LA
124  } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
125  edm::LogPrint("SiPixelLorentzAngleDB")
126  << " pixel endcap:"
127  << " side=" << tTopo->pxfSide(detidc.rawId()) << " disk=" << tTopo->pxfDisk(detidc.rawId())
128  << " blade=" << tTopo->pxfBlade(detidc.rawId()) << " panel=" << tTopo->pxfPanel(detidc.rawId())
129  << " module=" << tTopo->pxfModule(detidc.rawId()) << " rawId=" << detidc.rawId() << endl;
130 
131  //first individuals are put
132  for (Parameters::iterator it = ModuleParameters_.begin(); it != ModuleParameters_.end(); ++it) {
133  if (it->getParameter<unsigned int>("rawid") == detidc.rawId()) {
134  float lorentzangle = (float)it->getParameter<double>("angle");
135  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
136  edm::LogPrint("SiPixelLorentzAngleDB")
137  << " individual value=" << lorentzangle << " put into rawid=" << detid.rawId() << endl;
138  }
139  }
140 
141  //modules already put are automatically skipped
142  for (Parameters::iterator it = FPixParameters_.begin(); it != FPixParameters_.end(); ++it) {
143  if (it->getParameter<unsigned int>("side") == tTopo->pxfSide(detidc.rawId()) &&
144  it->getParameter<unsigned int>("disk") == tTopo->pxfDisk(detidc.rawId()) &&
145  it->getParameter<unsigned int>("HVgroup") ==
146  HVgroup(tTopo->pxfPanel(detidc.rawId()), tTopo->pxfModule(detidc.rawId()))) {
147  float lorentzangle = (float)it->getParameter<double>("angle");
148  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
149  }
150  }
151 
152  } else {
153  edm::LogError("SiPixelLorentzAngleDB")
154  << "[SiPixelLorentzAngleDB::analyze] detid is Pixel but neither bpix nor fpix" << std::endl;
155  }
156  }
157  }
158 
160  if (mydbservice.isAvailable()) {
161  try {
162  if (mydbservice->isNewTagRequest(recordName_)) {
163  mydbservice->createOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->beginOfTime(), recordName_);
164  } else {
165  mydbservice->appendOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->currentTime(), recordName_);
166  }
167  } catch (const cond::Exception& er) {
168  edm::LogError("SiPixelLorentzAngleDB") << er.what() << std::endl;
169  } catch (const std::exception& er) {
170  edm::LogError("SiPixelLorentzAngleDB") << "caught std::exception " << er.what() << std::endl;
171  } catch (...) {
172  edm::LogError("SiPixelLorentzAngleDB") << "Funny error" << std::endl;
173  }
174  } else {
175  edm::LogError("SiPixelLorentzAngleDB") << "Service is unavailable" << std::endl;
176  }
177 }
178 
179 unsigned int SiPixelLorentzAngleDB::HVgroup(unsigned int panel, unsigned int module) {
180  if (1 == panel && (1 == module || 2 == module)) {
181  return 1;
182  } else if (1 == panel && (3 == module || 4 == module)) {
183  return 2;
184  } else if (2 == panel && 1 == module) {
185  return 1;
186  } else if (2 == panel && (2 == module || 3 == module)) {
187  return 2;
188  } else {
189  edm::LogPrint("SiPixelLorentzAngleDB") << " *** error *** in SiPixelLorentzAngleDB::HVgroup(...), panel = " << panel
190  << ", module = " << module << endl;
191  return 0;
192  }
193 }
std::vector< std::pair< uint32_t, float > > detid_la
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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
unsigned int pxfBlade(const DetId &id) const
unsigned int pxfModule(const DetId &id) const
std::vector< edm::ParameterSet > Parameters
void analyze(const edm::Event &e, const edm::EventSetup &c) override
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
void appendOneIOV(const T &payload, cond::Time_t sinceTime, const std::string &recordName)
bool isNewTagRequest(const std::string &recordName)
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
unsigned int pxfDisk(const DetId &id) const
SiPixelLorentzAngleDB(const edm::ParameterSet &conf)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
unsigned int pxfPanel(const DetId &id) const
Log< level::Info, false > LogInfo
Definition: DetId.h:17
unsigned int pxfSide(const DetId &id) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_
unsigned int HVgroup(unsigned int panel, unsigned int module)
HLT enums.
bool isAvailable() const
Definition: Service.h:40
unsigned int pxbModule(const DetId &id) const
char const * what() const noexcept override
Definition: Exception.cc:107
~SiPixelLorentzAngleDB() override