Go to the documentation of this file.00001
00016 #include "RecoTBCalo/EcalTBRecProducers/interface/EcalTBWeightUncalibRecHitProducer.h"
00017 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00018 #include "DataFormats/EcalDigi/interface/EcalMGPASample.h"
00019 #include "DataFormats/Common/interface/Handle.h"
00020
00021 #include <iostream>
00022 #include <iomanip>
00023 #include <cmath>
00024
00025 #include "FWCore/Framework/interface/ESHandle.h"
00026
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028
00029 #include "CondFormats/EcalObjects/interface/EcalPedestals.h"
00030 #include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
00031
00032 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00033 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00034
00035 #include "CondFormats/EcalObjects/interface/EcalXtalGroupId.h"
00036 #include "CondFormats/EcalObjects/interface/EcalWeightXtalGroups.h"
00037 #include "CondFormats/DataRecord/interface/EcalWeightXtalGroupsRcd.h"
00038
00039 #include "CondFormats/EcalObjects/interface/EcalWeight.h"
00040 #include "CondFormats/EcalObjects/interface/EcalWeightSet.h"
00041 #include "CondFormats/EcalObjects/interface/EcalTBWeights.h"
00042 #include "CondFormats/DataRecord/interface/EcalTBWeightsRcd.h"
00043
00044 #include "CondFormats/EcalObjects/interface/EcalMGPAGainRatio.h"
00045 #include "CondFormats/EcalObjects/interface/EcalGainRatios.h"
00046 #include "CondFormats/DataRecord/interface/EcalGainRatiosRcd.h"
00047
00048 #include "CLHEP/Matrix/Matrix.h"
00049 #include "CLHEP/Matrix/SymMatrix.h"
00050 #include <vector>
00051
00052 #define DEBUG
00053 EcalTBWeightUncalibRecHitProducer::EcalTBWeightUncalibRecHitProducer(const edm::ParameterSet& ps) {
00054
00055 EBdigiCollection_ = ps.getParameter<edm::InputTag>("EBdigiCollection");
00056 EEdigiCollection_ = ps.getParameter<edm::InputTag>("EEdigiCollection");
00057 tdcRecInfoCollection_ = ps.getParameter<edm::InputTag>("tdcRecInfoCollection");
00058 EBhitCollection_ = ps.getParameter<std::string>("EBhitCollection");
00059 EEhitCollection_ = ps.getParameter<std::string>("EEhitCollection");
00060 nbTimeBin_ = ps.getParameter<int>("nbTimeBin");
00061 use2004OffsetConvention_ = ps.getUntrackedParameter< bool >("use2004OffsetConvention",false);
00062 produces< EBUncalibratedRecHitCollection >(EBhitCollection_);
00063 produces< EEUncalibratedRecHitCollection >(EEhitCollection_);
00064 }
00065
00066 EcalTBWeightUncalibRecHitProducer::~EcalTBWeightUncalibRecHitProducer() {
00067 }
00068
00069 void
00070 EcalTBWeightUncalibRecHitProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
00071
00072 using namespace edm;
00073
00074 Handle< EBDigiCollection > pEBDigis;
00075 const EBDigiCollection* EBdigis =0;
00076 if (EBdigiCollection_.label() != "" || EBdigiCollection_.instance() != "")
00077 {
00078
00079 evt.getByLabel( EBdigiCollection_, pEBDigis);
00080 if (!pEBDigis.isValid()) {
00081 edm::LogError("EcalUncalibRecHitError") << "Error! can't get the product " << EBdigiCollection_ ;
00082 } else {
00083 EBdigis = pEBDigis.product();
00084 #ifdef DEBUG
00085 LogDebug("EcalUncalibRecHitInfo") << "total # EBdigis: " << EBdigis->size() ;
00086 #endif
00087 }
00088 }
00089
00090 Handle< EEDigiCollection > pEEDigis;
00091 const EEDigiCollection* EEdigis =0;
00092
00093 if (EEdigiCollection_.label() != "" || EEdigiCollection_.instance() != "")
00094 {
00095
00096 evt.getByLabel( EEdigiCollection_, pEEDigis);
00097 if (!pEEDigis.isValid()) {
00098 edm::LogError("EcalUncalibRecHitError") << "Error! can't get the product " << EEdigiCollection_ ;
00099 } else {
00100 EEdigis = pEEDigis.product();
00101 #ifdef DEBUG
00102 LogDebug("EcalUncalibRecHitInfo") << "total # EEdigis: " << EEdigis->size() ;
00103
00104 #endif
00105 }
00106 }
00107
00108
00109 if (!EBdigis && !EEdigis)
00110 return;
00111
00112 Handle< EcalTBTDCRecInfo > pRecTDC;
00113 const EcalTBTDCRecInfo* recTDC =0;
00114
00115
00116 evt.getByLabel( tdcRecInfoCollection_, pRecTDC);
00117 if (pRecTDC.isValid()) {
00118 recTDC = pRecTDC.product();
00119 }
00120
00121
00122 edm::ESHandle<EcalWeightXtalGroups> pGrp;
00123 es.get<EcalWeightXtalGroupsRcd>().get(pGrp);
00124 const EcalWeightXtalGroups* grp = pGrp.product();
00125 if (!grp)
00126 return;
00127
00128 const EcalXtalGroupsMap& grpMap = grp->getMap();
00129
00130
00131
00132 edm::ESHandle<EcalGainRatios> pRatio;
00133 es.get<EcalGainRatiosRcd>().get(pRatio);
00134 const EcalGainRatioMap& gainMap = pRatio.product()->getMap();
00135
00136
00137
00138 #ifdef DEBUG
00139 LogDebug("EcalUncalibRecHitDebug") <<"Fetching EcalTBWeights from DB " ;
00140
00141 #endif
00142 edm::ESHandle<EcalTBWeights> pWgts;
00143 es.get<EcalTBWeightsRcd>().get(pWgts);
00144 const EcalTBWeights* wgts = pWgts.product();
00145
00146 if (!wgts)
00147 return;
00148
00149 #ifdef DEBUG
00150 LogDebug("EcalUncalibRecHitDebug") << "EcalTBWeightMap.size(): " << std::setprecision(3) << wgts->getMap().size() ;
00151
00152 #endif
00153
00154
00155 #ifdef DEBUG
00156 LogDebug("EcalUncalibRecHitDebug") << "fetching pedestals....";
00157
00158 #endif
00159 edm::ESHandle<EcalPedestals> pedHandle;
00160 es.get<EcalPedestalsRcd>().get( pedHandle );
00161 const EcalPedestalsMap& pedMap = pedHandle.product()->getMap();
00162 #ifdef DEBUG
00163 LogDebug("EcalUncalibRecHitDebug") << "done." ;
00164
00165 #endif
00166
00167
00168 std::auto_ptr< EBUncalibratedRecHitCollection > EBuncalibRechits( new EBUncalibratedRecHitCollection );
00169 std::auto_ptr< EEUncalibratedRecHitCollection > EEuncalibRechits( new EEUncalibratedRecHitCollection );
00170
00171 EcalPedestalsMapIterator pedIter;
00172
00173 EcalGainRatioMap::const_iterator gainIter;
00174
00175 EcalXtalGroupsMap::const_iterator git;
00176
00177 EcalTBWeights::EcalTBWeightMap::const_iterator wit;
00178
00179
00180 EcalTBWeights::EcalTDCId tdcid(int(nbTimeBin_/2)+1);
00181
00182 if (recTDC)
00183 if (recTDC->offset() == -999.)
00184 {
00185 edm::LogError("EcalUncalibRecHitError") << "TDC bin completely out of range. Returning" ;
00186 return;
00187 }
00188
00189 if (EBdigis)
00190 {
00191 for(unsigned int idig = 0; idig < EBdigis->size(); ++idig) {
00192
00193 EBDataFrame itdg = (*EBdigis)[idig];
00194
00195
00196
00197
00198
00199 #ifdef DEBUG
00200 LogDebug("EcalUncalibRecHitDebug") << "looking up pedestal for crystal: " << EBDetId(itdg.id()) ;
00201 #endif
00202 pedIter = pedMap.find(itdg.id().rawId());
00203 if( pedIter == pedMap.end() ) {
00204 edm::LogError("EcalUncalibRecHitError") << "error!! could not find pedestals for channel: " << EBDetId(itdg.id())
00205 << "\n no uncalib rechit will be made for this digi!"
00206 ;
00207 continue;
00208 }
00209 const EcalPedestals::Item& aped = (*pedIter);
00210 double pedVec[3];
00211 double pedRMSVec[3];
00212 pedVec[0]=aped.mean_x12;pedVec[1]=aped.mean_x6;pedVec[2]=aped.mean_x1;
00213 pedRMSVec[0]=aped.rms_x12;pedRMSVec[1]=aped.rms_x6;pedRMSVec[2]=aped.rms_x1;
00214
00215
00216 #ifdef DEBUG
00217 LogDebug("EcalUncalibRecHitDebug") << "looking up gainRatios for crystal: " << EBDetId(itdg.id()) ;
00218 #endif
00219 gainIter = gainMap.find(itdg.id().rawId());
00220 if( gainIter == gainMap.end() ) {
00221 edm::LogError("EcalUncalibRecHitError") << "error!! could not find gain ratios for channel: " << EBDetId(itdg.id())
00222 << "\n no uncalib rechit will be made for this digi!"
00223 ;
00224 continue;
00225 }
00226 const EcalMGPAGainRatio& aGain = (*gainIter);
00227 double gainRatios[3];
00228 gainRatios[0]=1.;gainRatios[1]=aGain.gain12Over6();gainRatios[2]=aGain.gain6Over1()*aGain.gain12Over6();
00229
00230
00231
00232
00233 git = grpMap.find( itdg.id().rawId() );
00234 if( git == grpMap.end() ) {
00235 edm::LogError("EcalUncalibRecHitError") << "No group id found for this crystal. something wrong with EcalWeightXtalGroups in your DB?"
00236 << "\n no uncalib rechit will be made for digi with id: " << EBDetId(itdg.id())
00237 ;
00238 continue;
00239 }
00240 const EcalXtalGroupId& gid = (*git);
00241
00242
00243
00244 double sampleGainRef = 1;
00245 int sampleSwitch = 999;
00246 for (int sample = 0; sample < itdg.size(); ++sample)
00247 {
00248 double gainSample = itdg.sample(sample).gainId();
00249 if(gainSample != sampleGainRef) {sampleGainRef = gainSample; sampleSwitch = sample;}
00250 }
00252
00253 if (recTDC)
00254 {
00255 int tdcBin=0;
00256 if (recTDC->offset() <= 0.)
00257 tdcBin = 1;
00258 if (recTDC->offset() >= 1.)
00259 tdcBin = nbTimeBin_;
00260 else
00261 tdcBin = int(recTDC->offset()*float(nbTimeBin_))+1;
00262
00263 if (tdcBin < 1 || tdcBin > nbTimeBin_ )
00264 {
00265 edm::LogError("EcalUncalibRecHitError") << "TDC bin out of range " << tdcBin << " offset " << recTDC->offset();
00266 continue;
00267 }
00268
00269
00270
00271
00272
00273 if (use2004OffsetConvention_ && sampleSwitch == 5)
00274 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00275 else if (!use2004OffsetConvention_ && sampleSwitch == 4)
00276 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00277 else
00278 tdcid=EcalTBWeights::EcalTDCId(tdcBin);
00279 }
00280
00281
00282 wit = wgts->getMap().find( std::make_pair(gid,tdcid) );
00283 if( wit == wgts->getMap().end() ) {
00284 edm::LogError("EcalUncalibRecHitError") << "No weights found for EcalGroupId: " << gid.id() << " and EcalTDCId: " << tdcid
00285 << "\n skipping digi with id: " << EBDetId(itdg.id())
00286 ;
00287 continue;
00288 }
00289 const EcalWeightSet& wset = wit->second;
00290
00291
00292
00293
00294 #ifdef DEBUG
00295 LogDebug("EcalUncalibRecHitDebug") << "accessing matrices of weights...";
00296 #endif
00297 const EcalWeightSet::EcalWeightMatrix& mat1 = wset.getWeightsBeforeGainSwitch();
00298 const EcalWeightSet::EcalWeightMatrix& mat2 = wset.getWeightsAfterGainSwitch();
00299
00300
00301 const EcalWeightSet::EcalWeightMatrix* weights[2];
00302 weights[0]=&mat1;
00303 weights[1]=&mat2;
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315 EcalUncalibratedRecHit aHit =
00316 EBalgo_.makeRecHit(itdg, pedVec, pedRMSVec, gainRatios, weights, testbeamEBShape);
00317
00318 EBuncalibRechits->push_back( aHit );
00319 #ifdef DEBUG
00320 if(aHit.amplitude()>0.) {
00321 LogDebug("EcalUncalibRecHitDebug") << "processed EBDataFrame with id: "
00322 << EBDetId(itdg.id()) << "\n"
00323 << "uncalib rechit amplitude: " << aHit.amplitude()
00324 ;
00325 }
00326 #endif
00327 }
00328 }
00329
00330 evt.put( EBuncalibRechits, EBhitCollection_ );
00331
00332
00333 if (EEdigis)
00334 {
00335 for(unsigned int idig = 0; idig < EEdigis->size(); ++idig) {
00336
00337 EEDataFrame itdg = (*EEdigis)[idig];
00338
00339
00340
00341
00342 #ifdef DEBUG
00343 LogDebug("EcalUncalibRecHitDebug") << "looking up pedestal for crystal: " << EEDetId(itdg.id()) ;
00344
00345 #endif
00346 pedIter = pedMap.find(itdg.id().rawId());
00347 if( pedIter == pedMap.end() ) {
00348 edm::LogError("EcalUncalibRecHitError") << "error!! could not find pedestals for channel: " << EEDetId(itdg.id())
00349 << "\n no uncalib rechit will be made for this digi!"
00350 ;
00351 continue;
00352 }
00353 const EcalPedestals::Item& aped = (*pedIter);
00354 double pedVec[3];
00355 double pedRMSVec[3];
00356 pedVec[0]=aped.mean_x12;pedVec[1]=aped.mean_x6;pedVec[2]=aped.mean_x1;
00357 pedRMSVec[0]=aped.rms_x12;pedRMSVec[1]=aped.rms_x6;pedRMSVec[2]=aped.rms_x1;
00358
00359
00360 #ifdef DEBUG
00361 LogDebug("EcalUncalibRecHitDebug") << "looking up gainRatios for crystal: " << EEDetId(itdg.id()) ;
00362
00363 #endif
00364 gainIter = gainMap.find(itdg.id().rawId());
00365 if( gainIter == gainMap.end() ) {
00366 edm::LogError("EcalUncalibRecHitError") << "error!! could not find gain ratios for channel: " << EEDetId(itdg.id())
00367 << "\n no uncalib rechit will be made for this digi!"
00368 ;
00369 continue;
00370 }
00371 const EcalMGPAGainRatio& aGain = (*gainIter);
00372 double gainRatios[3];
00373 gainRatios[0]=1.;gainRatios[1]=aGain.gain12Over6();gainRatios[2]=aGain.gain6Over1()*aGain.gain12Over6();
00374
00375
00376
00377
00378 git = grpMap.find( itdg.id().rawId() );
00379 if( git == grpMap.end() ) {
00380 edm::LogError("EcalUncalibRecHitError") << "No group id found for this crystal. something wrong with EcalWeightXtalGroups in your DB?"
00381 << "\n no uncalib rechit will be made for digi with id: " << EEDetId(itdg.id())
00382 ;
00383 continue;
00384 }
00385 const EcalXtalGroupId& gid = (*git);
00386
00387
00388
00389 double sampleGainRef = 1;
00390 int sampleSwitch = 999;
00391 for (int sample = 0; sample < itdg.size(); ++sample)
00392 {
00393 double gainSample = itdg.sample(sample).gainId();
00394 if(gainSample != sampleGainRef) {sampleGainRef = gainSample; sampleSwitch = sample;}
00395 }
00397
00398 if (recTDC)
00399 {
00400 int tdcBin=0;
00401 if (recTDC->offset() <= 0.)
00402 tdcBin = 1;
00403 if (recTDC->offset() >= 1.)
00404 tdcBin = nbTimeBin_;
00405 else
00406 tdcBin = int(recTDC->offset()*float(nbTimeBin_))+1;
00407
00408 if (tdcBin < 1 || tdcBin > nbTimeBin_ )
00409 {
00410 edm::LogError("EcalUncalibRecHitError") << "TDC bin out of range " << tdcBin << " offset " << recTDC->offset();
00411 continue;
00412 }
00413
00414
00415
00416
00417
00418 if (use2004OffsetConvention_ && sampleSwitch == 5)
00419 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00420 else if (!use2004OffsetConvention_ && sampleSwitch == 4)
00421 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00422 else
00423 tdcid=EcalTBWeights::EcalTDCId(tdcBin);
00424 }
00425
00426
00427 wit = wgts->getMap().find( std::make_pair(gid,tdcid) );
00428 if( wit == wgts->getMap().end() ) {
00429 edm::LogError("EcalUncalibRecHitError") << "No weights found for EcalGroupId: " << gid.id() << " and EcalTDCId: " << tdcid
00430 << "\n skipping digi with id: " << EEDetId(itdg.id())
00431 ;
00432 continue;
00433 }
00434 const EcalWeightSet& wset = wit->second;
00435
00436
00437
00438
00439 #ifdef DEBUG
00440 LogDebug("EcalUncalibRecHitDebug") << "accessing matrices of weights...";
00441
00442 #endif
00443 const EcalWeightSet::EcalWeightMatrix& mat1 = wset.getWeightsBeforeGainSwitch();
00444 const EcalWeightSet::EcalWeightMatrix& mat2 = wset.getWeightsAfterGainSwitch();
00445
00446
00447 const EcalWeightSet::EcalWeightMatrix* weights[2];
00448 weights[0]=&mat1;
00449 weights[1]=&mat2;
00450
00451
00452
00453
00454
00455
00456
00457
00458
00459
00460
00461 EcalUncalibratedRecHit aHit =
00462 EEalgo_.makeRecHit(itdg, pedVec, pedRMSVec, gainRatios, weights, testbeamEEShape);
00463
00464 EEuncalibRechits->push_back( aHit );
00465 #ifdef DEBUG
00466 if(aHit.amplitude()>0.) {
00467 LogDebug("EcalUncalibRecHitDebug") << "processed EEDataFrame with id: "
00468 << EEDetId(itdg.id()) << "\n"
00469 << "uncalib rechit amplitude: " << aHit.amplitude()
00470
00471
00472
00473
00474 ;
00475 }
00476 #endif
00477 }
00478 }
00479
00480 evt.put( EEuncalibRechits, EEhitCollection_ );
00481 }
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508