CMS 3D CMS Logo

SiPixelLorentzAngleDBLoader.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <iostream>
4 #include <fstream>
25 
26 class SiPixelLorentzAngleDBLoader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
27 public:
28  explicit SiPixelLorentzAngleDBLoader(const edm::ParameterSet& conf);
29 
30  ~SiPixelLorentzAngleDBLoader() override = default;
31  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
32 
33 private:
34  typedef std::vector<edm::ParameterSet> Parameters;
38  const bool useFile_;
45 
46  int HVgroup(int panel, int module);
47 
48  std::vector<std::pair<uint32_t, float> > detid_la;
49 };
50 
51 using namespace std;
52 using namespace edm;
53 
54 // Constructor
56  : tkGeomToken_(esConsumes()),
57  tkTopoToken_(esConsumes()),
58  recordName_(conf.getUntrackedParameter<std::string>("record", "SiPixelLorentzAngleRcd")),
59  useFile_(conf.getParameter<bool>("useFile")),
60  fileName_(conf.getParameter<string>("fileName")),
61  BPixParameters_(conf.getUntrackedParameter<Parameters>("BPixParameters")),
62  FPixParameters_(conf.getUntrackedParameter<Parameters>("FPixParameters")),
63  ModuleParameters_(conf.getUntrackedParameter<Parameters>("ModuleParameters")) {
65  static_cast<float>(conf.getUntrackedParameter<double>("bPixLorentzAnglePerTesla", -9999.));
67  static_cast<float>(conf.getUntrackedParameter<double>("fPixLorentzAnglePerTesla", -9999.));
69 }
70 
72  static constexpr int nModules_ = 4;
74 
75  // Retrieve tracker geometry from geometry
76  const TrackerGeometry* pDD = &es.getData(tkGeomToken_);
77  // Retrieve tracker topology from geometry
78  const TrackerTopology* tTopo = &es.getData(tkTopoToken_);
79 
80  for (auto& unit : pDD->detUnits()) {
81  if (auto pixelUnit = dynamic_cast<PixelGeomDetUnit const*>(unit)) {
82  const DetId detid = pixelUnit->geographicalId();
83  auto rawId = detid.rawId();
84  int found = 0;
85  int side = tTopo->side(detid); // 1:-z 2:+z for fpix, for bpix gives 0
86 
87  // fill bpix values for LA
88  if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
89  int layer = tTopo->pxbLayer(detid);
90  // Barrel ladder id 1-20,32,44.
91  int ladder = tTopo->pxbLadder(detid);
92  // Barrel Z-index=1,8
93  int module = tTopo->pxbModule(detid);
94  if (module < nModules_ + 1) {
95  side = 1;
96  } else {
97  side = 2;
98  }
99 
100  LogPrint("SiPixelLorentzAngleDBLoader") << " pixel barrel:"
101  << " layer=" << layer << " ladder=" << ladder << " module=" << module
102  << " rawId=" << rawId << " side=" << side;
103 
104  // use a commmon value (e.g. for MC)
105  if (bPixLorentzAnglePerTesla_ != -9999.) { // use common value for all
106  LogPrint("SiPixelLorentzAngleDBLoader")
107  << " LA=" << bPixLorentzAnglePerTesla_ << " common for all bpix" << endl;
108  if (!LorentzAngle.putLorentzAngle(detid.rawId(), bPixLorentzAnglePerTesla_)) {
109  LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
110  }
111  // use an external file
112  } else if (useFile_) {
113  LogPrint("SiPixelLorentzAngleDBLoader") << "method for reading file not implemented yet";
114  // use config file
115  } else {
116  // first individuals are put
117  for (auto& moduleParam : ModuleParameters_) {
118  if (moduleParam.getParameter<unsigned int>("rawid") == detid.rawId()) {
119  float lorentzangle = static_cast<float>(moduleParam.getParameter<double>("angle"));
120  if (!found) {
121  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
122  LogPrint("SiPixelLorentzAngleDBLoader")
123  << " >> LA=" << lorentzangle << " individual value " << detid.rawId();
124  found = 1;
125  } else {
126  LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
127  }
128  }
129  } // end on loop for ModuleParameters_
130 
131  //modules already put are automatically skipped
132  for (auto& bpixParam : BPixParameters_) {
133  if (bpixParam.exists("layer")) {
134  if (bpixParam.getParameter<int>("layer") != layer)
135  continue;
136  if (bpixParam.exists("ladder"))
137  if (bpixParam.getParameter<int>("ladder") != ladder)
138  continue;
139  if (bpixParam.exists("module"))
140  if (bpixParam.getParameter<int>("module") != module)
141  continue;
142  if (bpixParam.exists("side"))
143  if (bpixParam.getParameter<int>("side") != side)
144  continue;
145  if (!found) {
146  float lorentzangle = static_cast<float>(bpixParam.getParameter<double>("angle"));
147  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
148  LogPrint("SiPixelLorentzAngleDBLoader") << " >> LA=" << lorentzangle;
149  found = 2;
150  } else if (found == 1) {
151  LogPrint("SiPixelLorentzAngleDBLoader") << "The detid already given in ModuleParameters, skipping ...";
152  } else
153  LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
154  }
155  }
156  } // condition to read from config
157 
158  // fill fpix values for LA (for phase2 fpix & epix)
159  } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
160  // Convert to online
161  PixelEndcapName pen(detid, tTopo, true);
162  int disk = pen.diskName();
163  int blade = pen.bladeName();
164  int panel = pen.pannelName();
165  int ring = pen.ringName();
166 
167  LogPrint("SiPixelLorentzAngleDBLoader") << " pixel endcap:"
168  << " side=" << side << " disk=" << disk << " blade=" << blade
169  << " pannel=" << panel << " ring=" << ring << " rawId=" << rawId;
170 
171  // use a commmon value (e.g. for MC)
172  if (fPixLorentzAnglePerTesla_ != -9999.) { // use common value for all
173  LogPrint("SiPixelLorentzAngleDBLoader") << " LA =" << fPixLorentzAnglePerTesla_ << " common for all FPix";
174  if (!LorentzAngle.putLorentzAngle(detid.rawId(), fPixLorentzAnglePerTesla_)) {
175  LogError("SiPixelLorentzAngleDBLoader") << "detid already exists";
176  }
177 
178  } else if (useFile_) {
179  LogPrint("SiPixelLorentzAngleDBLoader") << "method for reading file not implemented yet";
180 
181  } else {
182  //first individuals are put
183  for (auto& parameter : ModuleParameters_) {
184  if (parameter.getParameter<unsigned int>("rawid") == detid.rawId()) {
185  float lorentzangle = static_cast<float>(parameter.getParameter<double>("angle"));
186  if (!found) {
187  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
188  LogPrint("SiPixelLorentzAngleDBLoader")
189  << " LA=" << lorentzangle << " individual value " << detid.rawId();
190  found = 1;
191  } else
192  LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
193  }
194  } // end loop on ModuleParameters_
195 
196  // modules already put are automatically skipped
197  for (auto& fpixParam : FPixParameters_) {
198  if (fpixParam.exists("side"))
199  if (fpixParam.getParameter<int>("side") != side)
200  continue;
201  if (fpixParam.exists("disk"))
202  if (fpixParam.getParameter<int>("disk") != disk)
203  continue;
204  if (fpixParam.exists("ring"))
205  if (fpixParam.getParameter<int>("ring") != ring)
206  continue;
207  if (fpixParam.exists("blade"))
208  if (fpixParam.getParameter<int>("blade") != blade)
209  continue;
210  if (fpixParam.exists("panel"))
211  if (fpixParam.getParameter<int>("panel") != panel)
212  continue;
213  if (fpixParam.exists("HVgroup"))
214  if (fpixParam.getParameter<int>("HVgroup") != HVgroup(panel, ring))
215  continue;
216  if (!found) {
217  float lorentzangle = static_cast<float>(fpixParam.getParameter<double>("angle"));
218  LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
219  LogPrint("SiPixelLorentzAngleDBLoader") << " >> LA=" << lorentzangle;
220  found = 2;
221  } else if (found == 1) {
222  LogPrint("SiPixelLorentzAngleDBLoader") << "The detid already given in ModuleParameters, skipping ...";
223  } else
224  LogError("SiPixelLorentzAngleDBLoader") << " ERROR!: detid already exists";
225  } // end loop on FPixParameters_
226  } // condition to read from config
227  } // end on being barrel or endcap
228  }
229  }
230 
232  if (mydbservice.isAvailable()) {
233  try {
234  if (mydbservice->isNewTagRequest(recordName_)) {
235  mydbservice->createOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->beginOfTime(), recordName_);
236  } else {
237  mydbservice->appendOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->currentTime(), recordName_);
238  }
239  } catch (const cond::Exception& er) {
240  LogPrint("SiPixelLorentzAngleDBLoader") << "SiPixelLorentzAngleDBLoader" << er.what();
241  } catch (const std::exception& er) {
242  LogPrint("SiPixelLorentzAngleDBLoader") << "SiPixelLorentzAngleDBLoader"
243  << "caught std::exception " << er.what();
244  }
245  } else {
246  LogPrint("SiPixelLorentzAngleDBLoader") << "Service is unavailable";
247  }
248 }
249 
251  if (1 == panel && (1 == module || 2 == module)) {
252  return 1;
253  } else if (1 == panel && (3 == module || 4 == module)) {
254  return 2;
255  } else if (2 == panel && 1 == module) {
256  return 1;
257  } else if (2 == panel && (2 == module || 3 == module)) {
258  return 2;
259  } else {
260  LogError("SiPixelLorentzAngleDBLoader")
261  << " *** error *** in SiPixelLorentzAngleDBLoader::HVgroup(...), panel = " << panel << ", module = " << module
262  << endl;
263  return 0;
264  }
265 }
266 
267 //define this as a plug-in
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
int bladeName() const
blade id
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
unsigned int pxbLadder(const DetId &id) const
Log< level::Error, false > LogError
unsigned int side(const DetId &id) const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
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)
int diskName() const
disk id
std::vector< std::pair< uint32_t, float > > detid_la
~SiPixelLorentzAngleDBLoader() override=default
SiPixelLorentzAngleDBLoader(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
Basic3DVector unit() const
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< edm::ParameterSet > Parameters
static const std::string kSharedResource
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
int pannelName() const
pannel id
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_