CMS 3D CMS Logo

EcalPFRecHitThresholdsMaker.cc
Go to the documentation of this file.
3 
5 
8 
10 
13 
47 
56 
59 
62 
71 
78 
81 
84 
87 
89 
91 
92 #include <vector>
93 
95  : m_timetype(iConfig.getParameter<std::string>("timetype")) {
96  std::string container;
99 
100  m_nsigma = iConfig.getParameter<double>("NSigma");
101 }
102 
104 
107  if (!dbOutput.isAvailable()) {
108  throw cms::Exception("PoolDBOutputService is not available");
109  }
110 
112  evtSetup.get<EcalPedestalsRcd>().get(handle1);
113  const EcalPedestals* ped_db = handle1.product();
114  std::cout << "ped pointer is: " << ped_db << std::endl;
115 
117  evtSetup.get<EcalADCToGeVConstantRcd>().get(handle2);
118  const EcalADCToGeVConstant* adc_db = handle2.product();
119  std::cout << "adc pointer is: " << adc_db << std::endl;
120 
122  evtSetup.get<EcalIntercalibConstantsRcd>().get(handle3);
123  const EcalIntercalibConstants* ical_db = handle3.product();
124  std::cout << "inter pointer is: " << ical_db << std::endl;
125 
127  evtSetup.get<EcalLaserDbRecord>().get(laser);
128 
130 
131  // const EcalIntercalibConstantMap& icalMap = ical_db->getMap();
132 
133  float adc_EB = float(adc_db->getEEValue());
134  float adc_EE = float(adc_db->getEBValue());
135 
136  //edm::Timestamp tsince;
137 
138  for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
139  if (iEta == 0)
140  continue;
141  for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
142  // make an EBDetId since we need EBDetId::rawId() to be used as the key for the pedestals
143  if (EBDetId::validDetId(iEta, iPhi)) {
144  EBDetId ebdetid(iEta, iPhi, EBDetId::ETAPHIMODE);
145  EcalPedestals::const_iterator it = ped_db->find(ebdetid.rawId());
146  EcalPedestals::Item aped = (*it);
147 
148  EcalIntercalibConstants::const_iterator itc = ical_db->find(ebdetid.rawId());
149  float calib = (*itc);
150 
151  // get laser coefficient
152  float lasercalib = 1.;
153  lasercalib = laser->getLaserCorrection(ebdetid, evt.time()); // TODO correct time
154 
155  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EB * lasercalib * m_nsigma;
156 
157  if (iPhi == 100)
158  std::cout << "Thresh(GeV)=" << thresh << std::endl;
159 
160  pfthresh->insert(std::make_pair(ebdetid.rawId(), thresh));
161  }
162  }
163  }
164 
165  for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
166  for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
167  // make an EEDetId since we need EEDetId::rawId() to be used as the key for the pedestals
168  if (EEDetId::validDetId(iX, iY, 1)) {
169  EEDetId eedetid(iX, iY, 1);
170 
171  EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
172  EcalPedestals::Item aped = (*it);
173 
174  EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
175  float calib = (*itc);
176 
177  // get laser coefficient
178  float lasercalib = 1.;
179  lasercalib = laser->getLaserCorrection(eedetid, evt.time()); // TODO correct time
180 
181  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
182  pfthresh->insert(std::make_pair(eedetid.rawId(), thresh));
183  }
184  if (EEDetId::validDetId(iX, iY, -1)) {
185  EEDetId eedetid(iX, iY, -1);
186 
187  EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
188  EcalPedestals::Item aped = (*it);
189 
190  EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
191  float calib = (*itc);
192 
193  // get laser coefficient
194  float lasercalib = 1.;
195  lasercalib = laser->getLaserCorrection(eedetid, evt.time()); // TODO correct time
196 
197  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
198  pfthresh->insert(std::make_pair(eedetid.rawId(), thresh));
199 
200  if (iX == 50)
201  std::cout << "Thresh(GeV)=" << thresh << std::endl;
202  }
203  }
204  }
205 
206  dbOutput->createNewIOV<const EcalPFRecHitThresholds>(
207  pfthresh, dbOutput->beginOfTime(), dbOutput->endOfTime(), "EcalPFRecHitThresholdsRcd");
208 
209  std::cout << "EcalPFRecHitThresholdsMaker wrote it " << std::endl;
210 }
T getParameter(std::string const &) const
static const int MIN_IPHI
Definition: EBDetId.h:135
EcalPFRecHitThresholdsMap EcalPFRecHitThresholds
EcalPFRecHitThresholdsMaker(const edm::ParameterSet &iConfig)
JetCorrectorParameters::Record record
Definition: classes.h:7
float getLaserCorrection(DetId const &xid, edm::Timestamp const &iTime) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
static const int IX_MIN
Definition: EEDetId.h:290
static const int IY_MIN
Definition: EEDetId.h:294
static bool validDetId(int i, int j)
check if a valid index combination
Definition: EBDetId.h:118
void analyze(const edm::Event &evt, const edm::EventSetup &evtSetup) override
float EcalPFRecHitThreshold
bool isAvailable() const
Definition: Service.h:40
static const int ETAPHIMODE
Definition: EBDetId.h:158
void createNewIOV(T *firstPayloadObj, cond::Time_t firstSinceTime, cond::Time_t firstTillTime, const std::string &recordName, bool withlogging=false)
void insert(std::pair< uint32_t, Item > const &a)
static const int IX_MAX
Definition: EEDetId.h:298
static const int MAX_IPHI
Definition: EBDetId.h:137
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
std::vector< Item >::const_iterator const_iterator
static const int MAX_IETA
Definition: EBDetId.h:136
T get() const
Definition: EventSetup.h:73
const_iterator find(uint32_t rawId) const
static const int IY_MAX
Definition: EEDetId.h:302
T const * product() const
Definition: ESHandle.h:86
edm::Timestamp time() const
Definition: EventBase.h:60