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 EcalTBWeightUncalibRecHitProducer::EcalTBWeightUncalibRecHitProducer(const edm::ParameterSet& ps) {
00053
00054 EBdigiCollection_ = ps.getParameter<std::string>("EBdigiCollection");
00055 digiProducer_ = ps.getParameter<std::string>("digiProducer");
00056 tdcRecInfoCollection_ = ps.getParameter<std::string>("tdcRecInfoCollection");
00057 tdcRecInfoProducer_ = ps.getParameter<std::string>("tdcRecInfoProducer");
00058 EBhitCollection_ = ps.getParameter<std::string>("EBhitCollection");
00059 nbTimeBin_ = ps.getParameter<int>("nbTimeBin");
00060 use2004OffsetConvention_ = ps.getUntrackedParameter< bool >("use2004OffsetConvention",false);
00061 produces< EBUncalibratedRecHitCollection >(EBhitCollection_);
00062 }
00063
00064 EcalTBWeightUncalibRecHitProducer::~EcalTBWeightUncalibRecHitProducer() {
00065 }
00066
00067 void
00068 EcalTBWeightUncalibRecHitProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
00069
00070 using namespace edm;
00071
00072 Handle< EBDigiCollection > pEBDigis;
00073 const EBDigiCollection* EBdigis =0;
00074
00075
00076 evt.getByLabel( digiProducer_, pEBDigis);
00077 if (!pEBDigis.isValid()) {
00078 edm::LogError("EcalUncalibRecHitError") << "Error! can't get the product " << EBdigiCollection_.c_str() ;
00079 } else {
00080 EBdigis = pEBDigis.product();
00081 #ifdef DEBUG
00082 LogDebug("EcalUncalibRecHitInfo") << "total # EBdigis: " << EBdigis->size() ;
00083 #endif
00084 }
00085
00086 if (!EBdigis)
00087 return;
00088
00089 Handle< EcalTBTDCRecInfo > pRecTDC;
00090 const EcalTBTDCRecInfo* recTDC =0;
00091
00092
00093 evt.getByLabel( tdcRecInfoProducer_, tdcRecInfoCollection_, pRecTDC);
00094 if (pRecTDC.isValid()) {
00095 recTDC = pRecTDC.product();
00096 }
00097
00098
00099 edm::ESHandle<EcalWeightXtalGroups> pGrp;
00100 es.get<EcalWeightXtalGroupsRcd>().get(pGrp);
00101 const EcalWeightXtalGroups* grp = pGrp.product();
00102 if (!grp)
00103 return;
00104
00105 const EcalXtalGroupsMap& grpMap = grp->getMap();
00106
00107
00108
00109 edm::ESHandle<EcalGainRatios> pRatio;
00110 es.get<EcalGainRatiosRcd>().get(pRatio);
00111 const EcalGainRatioMap& gainMap = pRatio.product()->getMap();
00112
00113
00114
00115 #ifdef DEBUG
00116 LogDebug("EcalUncalibRecHitDebug") <<"Fetching EcalTBWeights from DB " ;
00117 #endif
00118 edm::ESHandle<EcalTBWeights> pWgts;
00119 es.get<EcalTBWeightsRcd>().get(pWgts);
00120 const EcalTBWeights* wgts = pWgts.product();
00121
00122 if (!wgts)
00123 return;
00124
00125 #ifdef DEBUG
00126 LogDebug("EcalUncalibRecHitDebug") << "EcalTBWeightMap.size(): " << std::setprecision(3) << wgts->getMap().size() ;
00127 #endif
00128
00129
00130 #ifdef DEBUG
00131 LogDebug("EcalUncalibRecHitDebug") << "fetching pedestals....";
00132 #endif
00133 edm::ESHandle<EcalPedestals> pedHandle;
00134 es.get<EcalPedestalsRcd>().get( pedHandle );
00135 const EcalPedestalsMap& pedMap = pedHandle.product()->getMap();
00136 #ifdef DEBUG
00137 LogDebug("EcalUncalibRecHitDebug") << "done." ;
00138 #endif
00139
00140
00141 std::auto_ptr< EBUncalibratedRecHitCollection > EBuncalibRechits( new EBUncalibratedRecHitCollection );
00142
00143 EcalPedestalsMapIterator pedIter;
00144
00145 EcalGainRatioMap::const_iterator gainIter;
00146
00147 EcalXtalGroupsMap::const_iterator git;
00148
00149 EcalTBWeights::EcalTBWeightMap::const_iterator wit;
00150
00151
00152 EcalTBWeights::EcalTDCId tdcid(int(nbTimeBin_/2)+1);
00153
00154 if (recTDC)
00155 if (recTDC->offset() == -999.)
00156 {
00157 edm::LogError("EcalUncalibRecHitError") << "TDC bin completely out of range. Returning" ;
00158 return;
00159 }
00160
00161
00162 for(unsigned int idig = 0; idig < EBdigis->size(); ++idig) {
00163
00164 EBDataFrame itdg = (*EBdigis)[idig];
00165
00166
00167
00168
00169
00170 #ifdef DEBUG
00171 LogDebug("EcalUncalibRecHitDebug") << "looking up pedestal for crystal: " << EBDetId(itdg.id()) ;
00172 #endif
00173 pedIter = pedMap.find(itdg.id().rawId());
00174 if( pedIter == pedMap.end() ) {
00175 edm::LogError("EcalUncalibRecHitError") << "error!! could not find pedestals for channel: " << EBDetId(itdg.id())
00176 << "\n no uncalib rechit will be made for this digi!"
00177 ;
00178 continue;
00179 }
00180 const EcalPedestals::Item& aped = (*pedIter);
00181 double pedVec[3];
00182 pedVec[0]=aped.mean_x12;pedVec[1]=aped.mean_x6;pedVec[2]=aped.mean_x1;
00183
00184
00185 #ifdef DEBUG
00186 LogDebug("EcalUncalibRecHitDebug") << "looking up gainRatios for crystal: " << EBDetId(itdg.id()) ;
00187 #endif
00188 gainIter = gainMap.find(itdg.id().rawId());
00189 if( gainIter == gainMap.end() ) {
00190 edm::LogError("EcalUncalibRecHitError") << "error!! could not find gain ratios for channel: " << EBDetId(itdg.id())
00191 << "\n no uncalib rechit will be made for this digi!"
00192 ;
00193 continue;
00194 }
00195 const EcalMGPAGainRatio& aGain = (*gainIter);
00196 double gainRatios[3];
00197 gainRatios[0]=1.;gainRatios[1]=aGain.gain12Over6();gainRatios[2]=aGain.gain6Over1()*aGain.gain12Over6();
00198
00199
00200
00201
00202 git = grpMap.find( itdg.id().rawId() );
00203 if( git == grpMap.end() ) {
00204 edm::LogError("EcalUncalibRecHitError") << "No group id found for this crystal. something wrong with EcalWeightXtalGroups in your DB?"
00205 << "\n no uncalib rechit will be made for digi with id: " << EBDetId(itdg.id())
00206 ;
00207 continue;
00208 }
00209 const EcalXtalGroupId& gid = (*git);
00210
00211
00212
00213 double sampleGainRef = 1;
00214 int sampleSwitch = 999;
00215 for (int sample = 0; sample < itdg.size(); ++sample)
00216 {
00217 double gainSample = itdg.sample(sample).gainId();
00218 if(gainSample != sampleGainRef) {sampleGainRef = gainSample; sampleSwitch = sample;}
00219 }
00221
00222 if (recTDC)
00223 {
00224 int tdcBin=0;
00225 if (recTDC->offset() <= 0.)
00226 tdcBin = 1;
00227 if (recTDC->offset() >= 1.)
00228 tdcBin = nbTimeBin_;
00229 else
00230 tdcBin = int(recTDC->offset()*float(nbTimeBin_))+1;
00231
00232 if (tdcBin < 1 || tdcBin > nbTimeBin_ )
00233 {
00234 edm::LogError("EcalUncalibRecHitError") << "TDC bin out of range " << tdcBin << " offset " << recTDC->offset();
00235 continue;
00236 }
00237
00238
00239
00240
00241
00242 if (use2004OffsetConvention_ && sampleSwitch == 5)
00243 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00244 else if (!use2004OffsetConvention_ && sampleSwitch == 4)
00245 tdcid=EcalTBWeights::EcalTDCId(tdcBin+25);
00246 else
00247 tdcid=EcalTBWeights::EcalTDCId(tdcBin);
00248 }
00249
00250
00251 wit = wgts->getMap().find( std::make_pair(gid,tdcid) );
00252 if( wit == wgts->getMap().end() ) {
00253 edm::LogError("EcalUncalibRecHitError") << "No weights found for EcalGroupId: " << gid.id() << " and EcalTDCId: " << tdcid
00254 << "\n skipping digi with id: " << EBDetId(itdg.id())
00255 ;
00256 continue;
00257 }
00258 const EcalWeightSet& wset = wit->second;
00259
00260
00261
00262
00263 #ifdef DEBUG
00264 LogDebug("EcalUncalibRecHitDebug") << "accessing matrices of weights...";
00265 #endif
00266 const EcalWeightSet::EcalWeightMatrix& mat1 = wset.getWeightsBeforeGainSwitch();
00267 const EcalWeightSet::EcalWeightMatrix& mat2 = wset.getWeightsAfterGainSwitch();
00268 const EcalWeightSet::EcalChi2WeightMatrix& mat3 = wset.getChi2WeightsBeforeGainSwitch();
00269 const EcalWeightSet::EcalChi2WeightMatrix& mat4 = wset.getChi2WeightsAfterGainSwitch();
00270 const EcalWeightSet::EcalWeightMatrix* weights[2];
00271 weights[0]=&mat1;
00272 weights[1]=&mat2;
00273
00274
00275
00276
00277
00278
00279
00280 const EcalWeightSet::EcalChi2WeightMatrix* chi2mat[2];
00281 chi2mat[0]=&mat3;
00282 chi2mat[1]=&mat4;
00283
00284 EcalUncalibratedRecHit aHit =
00285 EBalgo_.makeRecHit(itdg, pedVec, gainRatios, weights, chi2mat);
00286 EBuncalibRechits->push_back( aHit );
00287 #ifdef DEBUG
00288 if(aHit.amplitude()>0.) {
00289 LogDebug("EcalUncalibRecHitDebug") << "processed EBDataFrame with id: "
00290 << EBDetId(itdg.id()) << "\n"
00291 << "uncalib rechit amplitude: " << aHit.amplitude()
00292 ;
00293 }
00294 #endif
00295 }
00296
00297 evt.put( EBuncalibRechits, EBhitCollection_ );
00298 }
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325