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