CMS 3D CMS Logo

EcalPFRecHitThresholdsMaker.cc
Go to the documentation of this file.
3 
5 
6 
9 
11 
14 
48 
57 
60 
63 
72 
79 
82 
85 
88 
90 
92 
93 #include <vector>
94 
96  m_timetype(iConfig.getParameter<std::string>("timetype"))
97 {
98 
99  std::string container;
102 
103  m_nsigma = iConfig.getParameter<double>("NSigma");
104 
105 }
106 
107 
109 {
110 
111 }
112 
114 
116  if ( !dbOutput.isAvailable() ) {
117  throw cms::Exception("PoolDBOutputService is not available");
118  }
119 
120 
121 
123  evtSetup.get<EcalPedestalsRcd>().get(handle1);
124  const EcalPedestals* ped_db = handle1.product();
125  std::cout << "ped pointer is: "<< ped_db<< std::endl;
126 
127 
129  evtSetup.get<EcalADCToGeVConstantRcd>().get(handle2);
130  const EcalADCToGeVConstant* adc_db = handle2.product();
131  std::cout << "adc pointer is: "<< adc_db<< std::endl;
132 
134  evtSetup.get<EcalIntercalibConstantsRcd>().get(handle3);
135  const EcalIntercalibConstants* ical_db = handle3.product();
136  std::cout << "inter pointer is: "<< ical_db<< std::endl;
137 
138 
140  evtSetup.get<EcalLaserDbRecord>().get(laser);
141 
143 
144  // const EcalIntercalibConstantMap& icalMap = ical_db->getMap();
145 
146 
147  float adc_EB= float(adc_db->getEEValue()) ;
148  float adc_EE=float(adc_db->getEBValue());
149 
150 
151  //edm::Timestamp tsince;
152 
153 
154 
155  for(int iEta=-EBDetId::MAX_IETA; iEta<=EBDetId::MAX_IETA ;++iEta) {
156  if(iEta==0) continue;
157  for(int iPhi=EBDetId::MIN_IPHI; iPhi<=EBDetId::MAX_IPHI; ++iPhi) {
158  // make an EBDetId since we need EBDetId::rawId() to be used as the key for the pedestals
159  if (EBDetId::validDetId(iEta,iPhi)) {
160  EBDetId ebdetid(iEta,iPhi,EBDetId::ETAPHIMODE);
161  EcalPedestals::const_iterator it =ped_db->find(ebdetid.rawId());
162  EcalPedestals::Item aped = (*it);
163 
164  EcalIntercalibConstants::const_iterator itc =ical_db->find(ebdetid.rawId());
165  float calib = (*itc);
166 
167  // get laser coefficient
168  float lasercalib = 1.;
169  lasercalib = laser->getLaserCorrection( ebdetid, evt.time()); // TODO correct time
170 
171  EcalPFRecHitThreshold thresh= aped.rms_x12 * calib * adc_EB * lasercalib * m_nsigma;
172 
173  if(iPhi==100) std::cout<<"Thresh(GeV)="<<thresh<<std::endl;
174 
175  pfthresh->insert(std::make_pair(ebdetid.rawId(),thresh));
176  }
177  }
178  }
179 
180  for(int iX=EEDetId::IX_MIN; iX<=EEDetId::IX_MAX ;++iX) {
181  for(int iY=EEDetId::IY_MIN; iY<=EEDetId::IY_MAX; ++iY) {
182  // make an EEDetId since we need EEDetId::rawId() to be used as the key for the pedestals
183  if (EEDetId::validDetId(iX,iY,1)) {
184  EEDetId eedetid(iX,iY,1);
185 
186  EcalPedestals::const_iterator it =ped_db->find(eedetid.rawId());
187  EcalPedestals::Item aped = (*it);
188 
189  EcalIntercalibConstants::const_iterator itc =ical_db->find(eedetid.rawId());
190  float calib = (*itc);
191 
192  // get laser coefficient
193  float lasercalib = 1.;
194  lasercalib = laser->getLaserCorrection( eedetid, evt.time()); // TODO correct time
195 
196  EcalPFRecHitThreshold thresh= aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
197  pfthresh->insert(std::make_pair(eedetid.rawId(),thresh));
198  }
199  if(EEDetId::validDetId(iX,iY,-1)) {
200 
201  EEDetId eedetid(iX,iY,-1);
202 
203  EcalPedestals::const_iterator it =ped_db->find(eedetid.rawId());
204  EcalPedestals::Item aped = (*it);
205 
206  EcalIntercalibConstants::const_iterator itc =ical_db->find(eedetid.rawId());
207  float calib = (*itc);
208 
209  // get laser coefficient
210  float lasercalib = 1.;
211  lasercalib = laser->getLaserCorrection( eedetid, evt.time()); // TODO correct time
212 
213  EcalPFRecHitThreshold thresh= aped.rms_x12 * calib * adc_EE * lasercalib * m_nsigma;
214  pfthresh->insert(std::make_pair(eedetid.rawId(),thresh));
215 
216  if(iX==50) std::cout<<"Thresh(GeV)="<<thresh<<std::endl;
217 
218 
219  }
220 
221  }
222  }
223 
224 
225 
226 
227  dbOutput->createNewIOV<const EcalPFRecHitThresholds>( pfthresh , dbOutput->beginOfTime(), dbOutput->endOfTime(),"EcalPFRecHitThresholdsRcd");
228 
229 
230 
231  std::cout<< "EcalPFRecHitThresholdsMaker wrote it " << std::endl;
232 }
233 
T getParameter(std::string const &) const
static const int MIN_IPHI
Definition: EBDetId.h:142
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:47
static const int IX_MIN
Definition: EEDetId.h:294
static const int IY_MIN
Definition: EEDetId.h:298
static bool validDetId(int i, int j)
check if a valid index combination
Definition: EBDetId.h:124
void analyze(const edm::Event &evt, const edm::EventSetup &evtSetup) override
float EcalPFRecHitThreshold
bool isAvailable() const
Definition: Service.h:46
static const int ETAPHIMODE
Definition: EBDetId.h:166
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:302
static const int MAX_IPHI
Definition: EBDetId.h:144
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:143
T get() const
Definition: EventSetup.h:63
const_iterator find(uint32_t rawId) const
static const int IY_MAX
Definition: EEDetId.h:306
T const * product() const
Definition: ESHandle.h:86
edm::Timestamp time() const
Definition: EventBase.h:61