CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripApvGainRescaler.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CondTools/SiStrip
4 // Class: SiStripApvGainRescaler
5 //
13 //
14 // Original Author: Marco Musich
15 // Created: Tue, 03 Oct 2017 12:57:34 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <iostream>
22 
23 // user include files
26 
29 
31 
36 
40 
41 //
42 // class declaration
43 //
44 
46 public:
48  ~SiStripApvGainRescaler() override;
49 
50  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
51 
52 private:
53  void beginJob() override;
54  void analyze(const edm::Event&, const edm::EventSetup&) override;
55  std::unique_ptr<SiStripApvGain> getNewObject(const std::map<std::pair<uint32_t, int>, float>& theMap);
56  void endJob() override;
57 
58  // ----------member data ---------------------------
60 
61  // take G2_old and G1_old from the regular gain handle
63  // take the additional G1_new from the Gain3Rcd (dirty trick)
65 };
66 
67 //
68 // constructors and destructor
69 //
71  : m_Record(iConfig.getParameter<std::string>("Record")), g1g2Token_(esConsumes()), g3Token_(esConsumes()) {
72  //now do what ever initialization is needed
73 }
74 
76 //
77 // member functions
78 //
79 
80 // ------------ method called for each event ------------
82  using namespace edm;
83 
84  const auto& g1g2 = iSetup.getData(g1g2Token_);
85  const auto& g3 = iSetup.getData(g3Token_);
86 
87  std::map<std::pair<uint32_t, int>, float> theMap;
88 
89  std::vector<uint32_t> detid;
90  g1g2.getDetIds(detid);
91  for (const auto& d : detid) {
92  SiStripApvGain::Range rangeG1_old = g1g2.getRange(d, 0);
93  SiStripApvGain::Range rangeG2_old = g1g2.getRange(d, 1);
94  SiStripApvGain::Range rangeG1_new = g3.getRange(d);
95 
96  int nAPV = 0;
97  for (int it = 0; it < rangeG1_old.second - rangeG1_old.first; it++) {
98  nAPV++;
99 
100  std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
101 
102  float G1_old = g1g2.getApvGain(it, rangeG1_old);
103  float G2_old = g1g2.getApvGain(it, rangeG2_old);
104  float G1G2_old = G1_old * G2_old;
105  float G1_new = g3.getApvGain(it, rangeG1_new);
106 
107  // this is based on G1_old*G2_old = G1_new * G2_new ==> G2_new = (G1_old*G2_old)/G1_new
108 
109  float NewGain = G1G2_old / G1_new;
110 
111  // DO NOT RESCALE APVs set to the default value
112  if (G2_old != 1.) {
113  theMap[index] = NewGain;
114  } else {
115  theMap[index] = 1.;
116  }
117 
118  } // loop over APVs
119  } // loop over DetIds
120 
121  std::unique_ptr<SiStripApvGain> theAPVGains = this->getNewObject(theMap);
122 
123  // write out the APVGains record
125 
126  if (poolDbService.isAvailable())
127  poolDbService->writeOneIOV(*theAPVGains, poolDbService->currentTime(), m_Record);
128  else
129  throw std::runtime_error("PoolDBService required.");
130 }
131 
132 // ------------ method called once each job just before starting event loop ------------
134 
135 // ------------ method called once each job just after ending the event loop ------------
137 
138 //********************************************************************************//
139 std::unique_ptr<SiStripApvGain> SiStripApvGainRescaler::getNewObject(
140  const std::map<std::pair<uint32_t, int>, float>& theMap) {
141  std::unique_ptr<SiStripApvGain> obj = std::make_unique<SiStripApvGain>();
142 
143  std::vector<float> theSiStripVector;
144  uint32_t PreviousDetId = 0;
145  for (const auto& element : theMap) {
146  uint32_t DetId = element.first.first;
147  if (DetId != PreviousDetId) {
148  if (!theSiStripVector.empty()) {
149  SiStripApvGain::Range range(theSiStripVector.begin(), theSiStripVector.end());
150  if (!obj->put(PreviousDetId, range))
151  printf("Bug to put detId = %i\n", PreviousDetId);
152  }
153  theSiStripVector.clear();
154  PreviousDetId = DetId;
155  }
156  theSiStripVector.push_back(element.second);
157 
158  edm::LogInfo("SiStripApvGainRescaler")
159  << " DetId: " << DetId << " APV: " << element.first.second << " Gain: " << element.second << std::endl;
160  }
161 
162  if (!theSiStripVector.empty()) {
163  SiStripApvGain::Range range(theSiStripVector.begin(), theSiStripVector.end());
164  if (!obj->put(PreviousDetId, range))
165  printf("Bug to put detId = %i\n", PreviousDetId);
166  }
167 
168  return obj;
169 }
170 
171 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
174 
175  desc.setComment(
176  " Utility class to rescale the values of SiStrip G2 by the ratio of G1_old/G1_new: this is useful in the case in "
177  "which a Gain2 payload needs to recycled after a G1 update to keep the G1*G2 product constant."
178  "PoolDBOutputService must be set up for 'SiStripApvGainRcd'.");
179 
180  desc.add<std::string>("Record", "SiStripApvGainRcd");
181  descriptions.add("rescaleGain2byGain1", desc);
182 }
183 
184 //define this as a plug-in
std::unique_ptr< SiStripApvGain > getNewObject(const std::map< std::pair< uint32_t, int >, float > &theMap)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SiStripApvGainRescaler(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const uint16_t range(const Frame &aFrame)
void analyze(const edm::Event &, const edm::EventSetup &) override
bool getData(T &iHolder) const
Definition: EventSetup.h:128
tuple d
Definition: ztail.py:151
void setComment(std::string const &value)
int iEvent
Definition: GenABIO.cc:224
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
std::pair< ContainerIterator, ContainerIterator > Range
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Info, false > LogInfo
Definition: DetId.h:17
~SiStripApvGainRescaler() override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::ESGetToken< SiStripApvGain, SiStripApvGain3Rcd > g3Token_
const edm::ESGetToken< SiStripGain, SiStripGainRcd > g1g2Token_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283