CMS 3D CMS Logo

EcalPFRecHitThresholdsMaker.cc
Go to the documentation of this file.
16 
17 #include <vector>
18 
20  : m_timetype(iConfig.getParameter<std::string>("timetype")),
21  ecalPedestalsToken_(esConsumes()),
22  ecalADCToGeVConstantToken_(esConsumes()),
23  ecalIntercalibConstantsToken_(esConsumes()),
24  ecalLaserDbServiceToken_(esConsumes()) {
25  m_nsigma = iConfig.getParameter<double>("NSigma");
26 }
27 
29 
32  if (!dbOutput.isAvailable()) {
33  throw cms::Exception("PoolDBOutputService is not available");
34  }
35 
36  const EcalPedestals* ped_db = &evtSetup.getData(ecalPedestalsToken_);
37  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "ped pointer is: " << ped_db << std::endl;
38 
39  const EcalADCToGeVConstant* adc_db = &evtSetup.getData(ecalADCToGeVConstantToken_);
40  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "adc pointer is: " << adc_db << std::endl;
41 
43  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "inter pointer is: " << ical_db << std::endl;
44 
45  const auto laser = evtSetup.getHandle(ecalLaserDbServiceToken_);
46 
47  EcalPFRecHitThresholds pfthresh;
48 
49  // const EcalIntercalibConstantMap& icalMap = ical_db->getMap();
50 
51  float adc_EB = float(adc_db->getEEValue());
52  float adc_EE = float(adc_db->getEBValue());
53 
54  //edm::Timestamp tsince;
55 
56  for (int iEta = -EBDetId::MAX_IETA; iEta <= EBDetId::MAX_IETA; ++iEta) {
57  if (iEta == 0)
58  continue;
59  for (int iPhi = EBDetId::MIN_IPHI; iPhi <= EBDetId::MAX_IPHI; ++iPhi) {
60  // make an EBDetId since we need EBDetId::rawId() to be used as the key for the pedestals
61  if (EBDetId::validDetId(iEta, iPhi)) {
62  EBDetId ebdetid(iEta, iPhi, EBDetId::ETAPHIMODE);
63  EcalPedestals::const_iterator it = ped_db->find(ebdetid.rawId());
64  EcalPedestals::Item aped = (*it);
65 
66  EcalIntercalibConstants::const_iterator itc = ical_db->find(ebdetid.rawId());
67  float calib = (*itc);
68 
69  // get laser coefficient
70  float lasercalib = 1.;
71  lasercalib = laser->getLaserCorrection(ebdetid, evt.time()); // TODO correct time
72 
73  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EB * lasercalib * m_nsigma;
74 
75  if (iPhi == 100)
76  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "Thresh(GeV)=" << thresh << std::endl;
77 
78  pfthresh.insert(std::make_pair(ebdetid.rawId(), thresh));
79  }
80  }
81  }
82 
83  for (int iX = EEDetId::IX_MIN; iX <= EEDetId::IX_MAX; ++iX) {
84  for (int iY = EEDetId::IY_MIN; iY <= EEDetId::IY_MAX; ++iY) {
85  // make an EEDetId since we need EEDetId::rawId() to be used as the key for the pedestals
86  if (EEDetId::validDetId(iX, iY, 1)) {
87  EEDetId eedetid(iX, iY, 1);
88 
89  EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
90  EcalPedestals::Item aped = (*it);
91 
92  EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
93  float calib = (*itc);
94 
95  // get laser coefficient
96  float lasercalib = 1.;
97  lasercalib = laser->getLaserCorrection(eedetid, evt.time()); // TODO correct time
98 
99  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
100  pfthresh.insert(std::make_pair(eedetid.rawId(), thresh));
101  }
102  if (EEDetId::validDetId(iX, iY, -1)) {
103  EEDetId eedetid(iX, iY, -1);
104 
105  EcalPedestals::const_iterator it = ped_db->find(eedetid.rawId());
106  EcalPedestals::Item aped = (*it);
107 
108  EcalIntercalibConstants::const_iterator itc = ical_db->find(eedetid.rawId());
109  float calib = (*itc);
110 
111  // get laser coefficient
112  float lasercalib = 1.;
113  lasercalib = laser->getLaserCorrection(eedetid, evt.time()); // TODO correct time
114 
115  EcalPFRecHitThreshold thresh = aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
116  pfthresh.insert(std::make_pair(eedetid.rawId(), thresh));
117 
118  if (iX == 50)
119  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "Thresh(GeV)=" << thresh << std::endl;
120  }
121  }
122  }
123 
124  dbOutput->createOneIOV<const EcalPFRecHitThresholds>(pfthresh, dbOutput->beginOfTime(), "EcalPFRecHitThresholdsRcd");
125 
126  edm::LogInfo("EcalPFRecHitThresholdsMaker") << "EcalPFRecHitThresholdsMaker wrote it " << std::endl;
127 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
static const int MIN_IPHI
Definition: EBDetId.h:135
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
EcalPFRecHitThresholdsMaker(const edm::ParameterSet &iConfig)
static const int IX_MIN
Definition: EEDetId.h:290
static const int IY_MIN
Definition: EEDetId.h:294
void createOneIOV(const T &payload, cond::Time_t firstSinceTime, const std::string &recordName)
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
edm::Timestamp time() const
Definition: EventBase.h:64
float EcalPFRecHitThreshold
edm::ESGetToken< EcalPedestals, EcalPedestalsRcd > ecalPedestalsToken_
edm::ESGetToken< EcalLaserDbService, EcalLaserDbRecord > ecalLaserDbServiceToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
static const int ETAPHIMODE
Definition: EBDetId.h:158
const_iterator find(uint32_t rawId) const
edm::ESGetToken< EcalIntercalibConstants, EcalIntercalibConstantsRcd > ecalIntercalibConstantsToken_
void insert(std::pair< uint32_t, Item > const &a)
static const int IX_MAX
Definition: EEDetId.h:298
Log< level::Info, false > LogInfo
static const int MAX_IPHI
Definition: EBDetId.h:137
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
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
edm::ESGetToken< EcalADCToGeVConstant, EcalADCToGeVConstantRcd > ecalADCToGeVConstantToken_
static const int IY_MAX
Definition: EEDetId.h:302
bool isAvailable() const
Definition: Service.h:40