00001
00002
00003
00004
00005
00006
00007 #include <map>
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "FWCore/Utilities/interface/EDMException.h"
00010 #include "FWCore/Framework/interface/ConstProductRegistry.h"
00011 #include "FWCore/ServiceRegistry/interface/Service.h"
00012 #include "FWCore/Framework/interface/ESHandle.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "DataFormats/Provenance/interface/Provenance.h"
00016 #include "DataFormats/Provenance/interface/BranchDescription.h"
00017
00018 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00019 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00020 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00021 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00022
00023
00024 #include "DataMixingHcalDigiWorker.h"
00025
00026
00027 using namespace std;
00028
00029 namespace edm
00030 {
00031
00032
00033
00034 DataMixingHcalDigiWorker::DataMixingHcalDigiWorker() { }
00035
00036
00037 DataMixingHcalDigiWorker::DataMixingHcalDigiWorker(const edm::ParameterSet& ps) :
00038 label_(ps.getParameter<std::string>("Label"))
00039
00040 {
00041
00042
00043
00044
00045
00046
00047
00048
00049 HBHEdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HBHEdigiCollectionSig");
00050 HOdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HOdigiCollectionSig");
00051 HFdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HFdigiCollectionSig");
00052 ZDCdigiCollectionSig_ = ps.getParameter<edm::InputTag>("ZDCdigiCollectionSig");
00053
00054 HBHEPileInputTag_ = ps.getParameter<edm::InputTag>("HBHEPileInputTag");
00055 HOPileInputTag_ = ps.getParameter<edm::InputTag>("HOPileInputTag");
00056 HFPileInputTag_ = ps.getParameter<edm::InputTag>("HFPileInputTag");
00057 ZDCPileInputTag_ = ps.getParameter<edm::InputTag>("ZDCPileInputTag");
00058
00059 DoZDC_ = false;
00060 if(ZDCPileInputTag_.label() != "") DoZDC_ = true;
00061
00062 HBHEDigiCollectionDM_ = ps.getParameter<std::string>("HBHEDigiCollectionDM");
00063 HODigiCollectionDM_ = ps.getParameter<std::string>("HODigiCollectionDM");
00064 HFDigiCollectionDM_ = ps.getParameter<std::string>("HFDigiCollectionDM");
00065 ZDCDigiCollectionDM_ = ps.getParameter<std::string>("ZDCDigiCollectionDM");
00066
00067
00068 }
00069
00070
00071 DataMixingHcalDigiWorker::~DataMixingHcalDigiWorker() {
00072 }
00073
00074 void DataMixingHcalDigiWorker::addHcalSignals(const edm::Event &e,const edm::EventSetup& ES) {
00075
00076
00077
00078 edm::ESHandle<HcalDbService> conditions;
00079 ES.get<HcalDbRecord>().get(conditions);
00080
00081
00082
00083
00084
00085 LogInfo("DataMixingHcalDigiWorker")<<"===============> adding MC signals for "<<e.id();
00086
00087
00088
00089 Handle< HBHEDigiCollection > pHBHEDigis;
00090
00091 const HBHEDigiCollection* HBHEDigis = 0;
00092
00093 if( e.getByLabel( HBHEdigiCollectionSig_, pHBHEDigis) ) {
00094 HBHEDigis = pHBHEDigis.product();
00095 LogDebug("DataMixingHcalDigiWorker") << "total # HBHE digis: " << HBHEDigis->size();
00096 }
00097
00098
00099
00100 if (HBHEDigis)
00101 {
00102
00103 for(HBHEDigiCollection::const_iterator it = HBHEDigis->begin();
00104 it != HBHEDigis->end(); ++it) {
00105
00106
00107
00108 HcalDetId cell = it->id();
00109
00110 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00111 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00112 HcalCoderDb coder (*channelCoder, *shape);
00113
00114 CaloSamples tool;
00115 coder.adc2fC((*it),tool);
00116
00117
00118
00119
00120
00121
00122
00123
00124 HBHEDigiStorage_.insert(HBHEDigiMap::value_type( ( it->id() ), tool ));
00125
00126 #ifdef DEBUG
00127 LogDebug("DataMixingHcalDigiWorker") << "processed HBHEDigi with rawId: "
00128 << it->id() << "\n"
00129 << " digi energy: " << it->energy();
00130 #endif
00131
00132 }
00133 }
00134
00135
00136
00137 Handle< HODigiCollection > pHODigis;
00138
00139 const HODigiCollection* HODigis = 0;
00140
00141 if( e.getByLabel( HOdigiCollectionSig_, pHODigis) ){
00142 HODigis = pHODigis.product();
00143 #ifdef DEBUG
00144 LogDebug("DataMixingHcalDigiWorker") << "total # HO digis: " << HODigis->size();
00145 #endif
00146 }
00147
00148
00149 if (HODigis)
00150 {
00151
00152 for(HODigiCollection::const_iterator it = HODigis->begin();
00153 it != HODigis->end(); ++it) {
00154
00155
00156 HcalDetId cell = it->id();
00157
00158 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00159 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00160 HcalCoderDb coder (*channelCoder, *shape);
00161
00162 CaloSamples tool;
00163 coder.adc2fC((*it),tool);
00164
00165 HODigiStorage_.insert(HODigiMap::value_type( ( it->id() ), tool ));
00166
00167 #ifdef DEBUG
00168 LogDebug("DataMixingHcalDigiWorker") << "processed HODigi with rawId: "
00169 << it->id() << "\n"
00170 << " digi energy: " << it->energy();
00171 #endif
00172
00173 }
00174 }
00175
00176
00177
00178 Handle< HFDigiCollection > pHFDigis;
00179
00180 const HFDigiCollection* HFDigis = 0;
00181
00182 if( e.getByLabel( HFdigiCollectionSig_, pHFDigis) ) {
00183 HFDigis = pHFDigis.product();
00184 #ifdef DEBUG
00185 LogDebug("DataMixingHcalDigiWorker") << "total # HF digis: " << HFDigis->size();
00186 #endif
00187 }
00188
00189
00190 if (HFDigis)
00191 {
00192
00193 for(HFDigiCollection::const_iterator it = HFDigis->begin();
00194 it != HFDigis->end(); ++it) {
00195
00196
00197 HcalDetId cell = it->id();
00198
00199 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00200 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00201 HcalCoderDb coder (*channelCoder, *shape);
00202
00203 CaloSamples tool;
00204 coder.adc2fC((*it),tool);
00205
00206 HFDigiStorage_.insert(HFDigiMap::value_type( ( it->id() ), tool ));
00207
00208 #ifdef DEBUG
00209 LogDebug("DataMixingHcalDigiWorker") << "processed HFDigi with rawId: "
00210 << it->id() << "\n"
00211 << " digi energy: " << it->energy();
00212 #endif
00213
00214 }
00215 }
00216
00217
00218
00219 if(DoZDC_){
00220
00221 Handle< ZDCDigiCollection > pZDCDigis;
00222
00223 const ZDCDigiCollection* ZDCDigis = 0;
00224
00225 if( e.getByLabel( ZDCdigiCollectionSig_, pZDCDigis) ) {
00226 ZDCDigis = pZDCDigis.product();
00227 #ifdef DEBUG
00228 LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
00229 #endif
00230 }
00231
00232
00233 if (ZDCDigis)
00234 {
00235
00236 for(ZDCDigiCollection::const_iterator it = ZDCDigis->begin();
00237 it != ZDCDigis->end(); ++it) {
00238
00239
00240 HcalDetId cell = it->id();
00241
00242 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00243 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00244 HcalCoderDb coder (*channelCoder, *shape);
00245
00246 CaloSamples tool;
00247 coder.adc2fC((*it),tool);
00248
00249 ZDCDigiStorage_.insert(ZDCDigiMap::value_type( ( it->id() ), tool ));
00250
00251 #ifdef DEBUG
00252 LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with rawId: "
00253 << it->id() << "\n"
00254 << " digi energy: " << it->energy();
00255 #endif
00256
00257 }
00258 }
00259 }
00260
00261 }
00262
00263 void DataMixingHcalDigiWorker::addHcalPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr,const edm::EventSetup& ES) {
00264
00265 LogDebug("DataMixingHcalDigiWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr;
00266
00267
00268 edm::ESHandle<HcalDbService> conditions;
00269 ES.get<HcalDbRecord>().get(conditions);
00270
00271
00272
00273
00274
00275 boost::shared_ptr<Wrapper<HBHEDigiCollection> const> HBHEDigisPTR =
00276 getProductByTag<HBHEDigiCollection>(*ep, HBHEPileInputTag_ );
00277
00278 if(HBHEDigisPTR ) {
00279
00280 const HBHEDigiCollection* HBHEDigis = const_cast< HBHEDigiCollection * >(HBHEDigisPTR->product());
00281
00282 LogInfo("DataMixingHcalDigiWorker") << "total # HBHE digis: " << HBHEDigis->size();
00283
00284
00285 for(HBHEDigiCollection::const_iterator it = HBHEDigis->begin();
00286 it != HBHEDigis->end(); ++it) {
00287
00288
00289 HcalDetId cell = it->id();
00290
00291 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00292 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00293 HcalCoderDb coder (*channelCoder, *shape);
00294
00295 CaloSamples tool;
00296 coder.adc2fC((*it),tool);
00297
00298 HBHEDigiStorage_.insert(HBHEDigiMap::value_type( (it->id()), tool ));
00299
00300 #ifdef DEBUG
00301 LogDebug("DataMixingHcalDigiWorker") << "processed HBHEDigi with rawId: "
00302 << it->id() << "\n"
00303 << " digi energy: " << it->energy();
00304 #endif
00305 }
00306 }
00307
00308
00309 boost::shared_ptr<Wrapper<HODigiCollection> const> HODigisPTR =
00310 getProductByTag<HODigiCollection>(*ep, HOPileInputTag_ );
00311
00312 if(HODigisPTR ) {
00313
00314 const HODigiCollection* HODigis = const_cast< HODigiCollection * >(HODigisPTR->product());
00315
00316 LogDebug("DataMixingHcalDigiWorker") << "total # HO digis: " << HODigis->size();
00317
00318
00319 for(HODigiCollection::const_iterator it = HODigis->begin();
00320 it != HODigis->end(); ++it) {
00321
00322
00323 HcalDetId cell = it->id();
00324
00325 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00326 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00327 HcalCoderDb coder (*channelCoder, *shape);
00328
00329 CaloSamples tool;
00330 coder.adc2fC((*it),tool);
00331
00332 HODigiStorage_.insert(HODigiMap::value_type( (it->id()), tool ));
00333
00334 #ifdef DEBUG
00335 LogDebug("DataMixingHcalDigiWorker") << "processed HODigi with rawId: "
00336 << it->id() << "\n"
00337 << " digi energy: " << it->energy();
00338 #endif
00339 }
00340 }
00341
00342
00343
00344
00345 boost::shared_ptr<Wrapper<HFDigiCollection> const> HFDigisPTR =
00346 getProductByTag<HFDigiCollection>(*ep, HFPileInputTag_ );
00347
00348 if(HFDigisPTR ) {
00349
00350 const HFDigiCollection* HFDigis = const_cast< HFDigiCollection * >(HFDigisPTR->product());
00351
00352 LogDebug("DataMixingHcalDigiWorker") << "total # HF digis: " << HFDigis->size();
00353
00354
00355 for(HFDigiCollection::const_iterator it = HFDigis->begin();
00356 it != HFDigis->end(); ++it) {
00357
00358
00359 HcalDetId cell = it->id();
00360
00361 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00362 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00363 HcalCoderDb coder (*channelCoder, *shape);
00364
00365 CaloSamples tool;
00366 coder.adc2fC((*it),tool);
00367
00368 HFDigiStorage_.insert(HFDigiMap::value_type( (it->id()), tool ));
00369
00370 #ifdef DEBUG
00371 LogDebug("DataMixingHcalDigiWorker") << "processed HFDigi with rawId: "
00372 << it->id() << "\n"
00373 << " digi energy: " << it->energy();
00374 #endif
00375 }
00376 }
00377
00378
00379
00380
00381 if(DoZDC_) {
00382
00383
00384 boost::shared_ptr<Wrapper<ZDCDigiCollection> const> ZDCDigisPTR =
00385 getProductByTag<ZDCDigiCollection>(*ep, ZDCPileInputTag_ );
00386
00387 if(ZDCDigisPTR ) {
00388
00389 const ZDCDigiCollection* ZDCDigis = const_cast< ZDCDigiCollection * >(ZDCDigisPTR->product());
00390
00391 LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
00392
00393
00394 for(ZDCDigiCollection::const_iterator it = ZDCDigis->begin();
00395 it != ZDCDigis->end(); ++it) {
00396
00397
00398 HcalDetId cell = it->id();
00399
00400 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00401 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00402 HcalCoderDb coder (*channelCoder, *shape);
00403
00404 CaloSamples tool;
00405 coder.adc2fC((*it),tool);
00406
00407 ZDCDigiStorage_.insert(ZDCDigiMap::value_type( (it->id()), tool ));
00408
00409 #ifdef DEBUG
00410 LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with rawId: "
00411 << it->id() << "\n"
00412 << " digi energy: " << it->energy();
00413 #endif
00414 }
00415 }
00416 }
00417
00418
00419 }
00420
00421 void DataMixingHcalDigiWorker::putHcal(edm::Event &e,const edm::EventSetup& ES) {
00422
00423
00424 std::auto_ptr< HBHEDigiCollection > HBHEdigis( new HBHEDigiCollection );
00425 std::auto_ptr< HODigiCollection > HOdigis( new HODigiCollection );
00426 std::auto_ptr< HFDigiCollection > HFdigis( new HFDigiCollection );
00427 std::auto_ptr< ZDCDigiCollection > ZDCdigis( new ZDCDigiCollection );
00428
00429
00430 edm::ESHandle<HcalDbService> conditions;
00431 ES.get<HcalDbRecord>().get(conditions);
00432
00433
00434 DetId formerID = 0;
00435 DetId currentID;
00436
00437 CaloSamples HB_old;
00438
00439 double fC_new;
00440 double fC_old;
00441 double fC_sum;
00442
00443
00444
00445 HBHEDigiMap::const_iterator iHBchk;
00446
00447 for(HBHEDigiMap::const_iterator iHB = HBHEDigiStorage_.begin();
00448 iHB != HBHEDigiStorage_.end(); ++iHB) {
00449
00450 currentID = iHB->first;
00451
00452 if (currentID == formerID) {
00453
00454
00455
00456 unsigned int sizenew = (iHB->second).size();
00457 unsigned int sizeold = HB_old.size();
00458
00459 bool usenew = false;
00460
00461 if(sizenew > sizeold) usenew = true;
00462
00463 unsigned int max_samp = std::max(sizenew, sizeold);
00464
00465 CaloSamples HB_bigger(currentID,max_samp);
00466
00467
00468
00469
00470
00471
00472 for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00473 if(isamp < sizenew) {
00474 fC_new = (iHB->second)[isamp];
00475 }
00476 else { fC_new = 0;}
00477
00478 if(isamp < sizeold) {
00479 fC_old = HB_old[isamp];
00480 }
00481 else { fC_old = 0;}
00482
00483
00484 fC_sum = fC_new + fC_old;
00485
00486
00487
00488
00489
00490 if(usenew) {HB_bigger[isamp] = fC_sum; }
00491 else { HB_old[isamp] = fC_sum; }
00492
00493 }
00494 if(usenew) HB_old = HB_bigger;
00495
00496 }
00497 else {
00498 if(formerID>0) {
00499
00500 HBHEdigis->push_back(HBHEDataFrame(formerID));
00501
00502
00503
00504 HcalDetId cell = HB_old.id();
00505 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00506 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00507 HcalCoderDb coder (*channelCoder, *shape);
00508
00509 unsigned int sizeold = HB_old.size();
00510 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00511 coder.fC2adc(HB_old,(HBHEdigis->back()), 0 );
00512 }
00513 }
00514
00515 formerID = currentID;
00516 HB_old = iHB->second;
00517
00518 }
00519
00520 iHBchk = iHB;
00521 if((++iHBchk) == HBHEDigiStorage_.end()) {
00522
00523
00524 HBHEdigis->push_back(HBHEDataFrame(currentID));
00525
00526
00527
00528 HcalDetId cell = (iHB->second).id();
00529 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00530 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00531 HcalCoderDb coder (*channelCoder, *shape);
00532
00533 unsigned int sizenew = (iHB->second).size();
00534 for(unsigned int isamp = 0; isamp<sizenew; isamp++) {
00535 coder.fC2adc(HB_old,(HBHEdigis->back()), 0 );
00536 }
00537 }
00538 }
00539
00540
00541
00542
00543
00544 formerID = 0;
00545 CaloSamples HO_old;
00546
00547 HODigiMap::const_iterator iHOchk;
00548
00549 for(HODigiMap::const_iterator iHO = HODigiStorage_.begin();
00550 iHO != HODigiStorage_.end(); ++iHO) {
00551
00552 currentID = iHO->first;
00553
00554 if (currentID == formerID) {
00555
00556
00557 unsigned int sizenew = (iHO->second).size();
00558 unsigned int sizeold = HO_old.size();
00559
00560 unsigned int max_samp = std::max(sizenew, sizeold);
00561
00562 CaloSamples HO_bigger(currentID,max_samp);
00563
00564 bool usenew = false;
00565
00566 if(sizenew > sizeold) usenew = true;
00567
00568
00569
00570
00571 for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00572 if(isamp < sizenew) {
00573 fC_new = (iHO->second)[isamp];
00574 }
00575 else { fC_new = 0;}
00576
00577 if(isamp < sizeold) {
00578 fC_old = HO_old[isamp];
00579 }
00580 else { fC_old = 0;}
00581
00582
00583 fC_sum = fC_new + fC_old;
00584
00585 if(usenew) {HO_bigger[isamp] = fC_sum; }
00586 else { HO_old[isamp] = fC_sum; }
00587
00588 }
00589 if(usenew) HO_old = HO_bigger;
00590
00591 }
00592 else {
00593 if(formerID>0) {
00594
00595 HOdigis->push_back(HODataFrame(formerID));
00596
00597
00598
00599 HcalDetId cell = HO_old.id();
00600 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00601 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00602 HcalCoderDb coder (*channelCoder, *shape);
00603
00604 unsigned int sizeold = HO_old.size();
00605 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00606 coder.fC2adc(HO_old,(HOdigis->back()), 0 );
00607 }
00608 }
00609
00610 formerID = currentID;
00611 HO_old = iHO->second;
00612 }
00613
00614 iHOchk = iHO;
00615 if((++iHOchk) == HODigiStorage_.end()) {
00616
00617 HOdigis->push_back(HODataFrame(currentID));
00618
00619
00620
00621 HcalDetId cell = (iHO->second).id();
00622 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00623 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00624 HcalCoderDb coder (*channelCoder, *shape);
00625
00626 unsigned int sizeold = (iHO->second).size();
00627 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00628 coder.fC2adc(HO_old,(HOdigis->back()), 0 );
00629 }
00630
00631 }
00632 }
00633
00634
00635
00636
00637 formerID = 0;
00638 CaloSamples HF_old;
00639
00640 HFDigiMap::const_iterator iHFchk;
00641
00642 for(HFDigiMap::const_iterator iHF = HFDigiStorage_.begin();
00643 iHF != HFDigiStorage_.end(); ++iHF) {
00644
00645 currentID = iHF->first;
00646
00647 if (currentID == formerID) {
00648
00649
00650 unsigned int sizenew = (iHF->second).size();
00651 unsigned int sizeold = HF_old.size();
00652
00653 unsigned int max_samp = std::max(sizenew, sizeold);
00654
00655 CaloSamples HF_bigger(currentID,max_samp);
00656
00657 bool usenew = false;
00658
00659 if(sizenew > sizeold) usenew = true;
00660
00661
00662
00663
00664 for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00665 if(isamp < sizenew) {
00666 fC_new = (iHF->second)[isamp];
00667 }
00668 else { fC_new = 0;}
00669
00670 if(isamp < sizeold) {
00671 fC_old = HF_old[isamp];
00672 }
00673 else { fC_old = 0;}
00674
00675
00676 fC_sum = fC_new + fC_old;
00677
00678 if(usenew) {HF_bigger[isamp] = fC_sum; }
00679 else { HF_old[isamp] = fC_sum; }
00680
00681 }
00682 if(usenew) HF_old = HF_bigger;
00683
00684 }
00685 else {
00686 if(formerID>0) {
00687
00688 HFdigis->push_back(HFDataFrame(formerID));
00689
00690
00691
00692 HcalDetId cell = HF_old.id();
00693 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00694 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00695 HcalCoderDb coder (*channelCoder, *shape);
00696
00697 unsigned int sizeold = HF_old.size();
00698 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00699 coder.fC2adc(HF_old,(HFdigis->back()), 0 );
00700 }
00701 }
00702
00703 formerID = currentID;
00704 HF_old = iHF->second;
00705 }
00706
00707 iHFchk = iHF;
00708 if((++iHFchk) == HFDigiStorage_.end()) {
00709
00710 HFdigis->push_back(HFDataFrame(currentID));
00711
00712
00713
00714 HcalDetId cell = (iHF->second).id();
00715 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00716 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00717 HcalCoderDb coder (*channelCoder, *shape);
00718
00719 unsigned int sizeold = (iHF->second).size();
00720 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00721 coder.fC2adc(HF_old,(HFdigis->back()), 0 );
00722 }
00723
00724 }
00725 }
00726
00727
00728
00729
00730
00731 formerID = 0;
00732 CaloSamples ZDC_old;
00733
00734 ZDCDigiMap::const_iterator iZDCchk;
00735
00736 for(ZDCDigiMap::const_iterator iZDC = ZDCDigiStorage_.begin();
00737 iZDC != ZDCDigiStorage_.end(); ++iZDC) {
00738
00739 currentID = iZDC->first;
00740
00741 if (currentID == formerID) {
00742
00743
00744 unsigned int sizenew = (iZDC->second).size();
00745 unsigned int sizeold = ZDC_old.size();
00746
00747 unsigned int max_samp = std::max(sizenew, sizeold);
00748
00749 CaloSamples ZDC_bigger(currentID,max_samp);
00750
00751 bool usenew = false;
00752
00753 if(sizenew > sizeold) usenew = true;
00754
00755
00756
00757
00758 for(unsigned int isamp = 0; isamp<max_samp; isamp++) {
00759 if(isamp < sizenew) {
00760 fC_new = (iZDC->second)[isamp];
00761 }
00762 else { fC_new = 0;}
00763
00764 if(isamp < sizeold) {
00765 fC_old = ZDC_old[isamp];
00766 }
00767 else { fC_old = 0;}
00768
00769
00770 fC_sum = fC_new + fC_old;
00771
00772 if(usenew) {ZDC_bigger[isamp] = fC_sum; }
00773 else { ZDC_old[isamp] = fC_sum; }
00774
00775 }
00776 if(usenew) ZDC_old = ZDC_bigger;
00777
00778 }
00779 else {
00780 if(formerID>0) {
00781
00782 ZDCdigis->push_back(ZDCDataFrame(formerID));
00783
00784
00785
00786 HcalDetId cell = ZDC_old.id();
00787 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00788 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00789 HcalCoderDb coder (*channelCoder, *shape);
00790
00791 unsigned int sizeold = ZDC_old.size();
00792 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00793 coder.fC2adc(ZDC_old,(ZDCdigis->back()), 0 );
00794 }
00795 }
00796
00797 formerID = currentID;
00798 ZDC_old = iZDC->second;
00799 }
00800
00801 iZDCchk = iZDC;
00802 if((++iZDCchk) == ZDCDigiStorage_.end()) {
00803
00804 ZDCdigis->push_back(ZDCDataFrame(currentID));
00805
00806
00807
00808 HcalDetId cell = (iZDC->second).id();
00809 const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00810 const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
00811 HcalCoderDb coder (*channelCoder, *shape);
00812
00813 unsigned int sizeold = (iZDC->second).size();
00814 for(unsigned int isamp = 0; isamp<sizeold; isamp++) {
00815 coder.fC2adc(ZDC_old,(ZDCdigis->back()), 0 );
00816 }
00817
00818 }
00819 }
00820
00821
00822
00823
00824
00825
00826 LogInfo("DataMixingHcalDigiWorker") << "total # HBHE Merged digis: " << HBHEdigis->size() ;
00827 LogInfo("DataMixingHcalDigiWorker") << "total # HO Merged digis: " << HOdigis->size() ;
00828 LogInfo("DataMixingHcalDigiWorker") << "total # HF Merged digis: " << HFdigis->size() ;
00829 LogInfo("DataMixingHcalDigiWorker") << "total # ZDC Merged digis: " << ZDCdigis->size() ;
00830
00831 e.put( HBHEdigis, HBHEDigiCollectionDM_ );
00832 e.put( HOdigis, HODigiCollectionDM_ );
00833 e.put( HFdigis, HFDigiCollectionDM_ );
00834 e.put( ZDCdigis, ZDCDigiCollectionDM_ );
00835
00836
00837 HBHEDigiStorage_.clear();
00838 HODigiStorage_.clear();
00839 HFDigiStorage_.clear();
00840 ZDCDigiStorage_.clear();
00841
00842 }
00843
00844 }