CMS 3D CMS Logo

DataMixingHcalWorker.cc

Go to the documentation of this file.
00001 // File: DataMixingHcalWorker.cc
00002 // Description:  see DataMixingHcalWorker.h
00003 // Author:  Mike Hildreth, University of Notre Dame
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 "DataFormats/Common/interface/Handle.h"
00013 #include "DataFormats/Provenance/interface/Provenance.h"
00014 #include "DataFormats/Provenance/interface/BranchDescription.h"
00015 //
00016 //
00017 #include "DataMixingHcalWorker.h"
00018 
00019 
00020 using namespace std;
00021 
00022 namespace edm
00023 {
00024 
00025   // Virtual constructor
00026 
00027   DataMixingHcalWorker::DataMixingHcalWorker() {sel_=0;} 
00028 
00029   // Constructor 
00030   DataMixingHcalWorker::DataMixingHcalWorker(const edm::ParameterSet& ps) : 
00031                                                             label_(ps.getParameter<std::string>("Label"))
00032 
00033   {                                                         
00034 
00035     // get the subdetector names
00036     //    this->getSubdetectorNames();  //something like this may be useful to check what we are supposed to do...
00037 
00038     // create input selector
00039     if (label_.size()>0){
00040       sel_=new Selector( ModuleLabelSelector(label_));
00041     }
00042     else {
00043       sel_=new Selector( MatchAllSelector());
00044     }
00045 
00046     // declare the products to produce
00047 
00048     // Hcal 
00049 
00050     HBHErechitCollection_  = ps.getParameter<edm::InputTag>("HBHEProducer");
00051     HOrechitCollection_    = ps.getParameter<edm::InputTag>("HOProducer");
00052     HFrechitCollection_    = ps.getParameter<edm::InputTag>("HFProducer");
00053     ZDCrechitCollection_   = ps.getParameter<edm::InputTag>("ZDCrechitCollection");
00054 
00055     HBHERecHitCollectionDM_ = ps.getParameter<std::string>("HBHERecHitCollectionDM");
00056     HORecHitCollectionDM_   = ps.getParameter<std::string>("HORecHitCollectionDM");
00057     HFRecHitCollectionDM_   = ps.getParameter<std::string>("HFRecHitCollectionDM");
00058     ZDCRecHitCollectionDM_  = ps.getParameter<std::string>("ZDCRecHitCollectionDM");
00059 
00060 
00061   }
00062                
00063   // Virtual destructor needed.
00064   DataMixingHcalWorker::~DataMixingHcalWorker() { 
00065     delete sel_;
00066     sel_=0;
00067   }  
00068 
00069   void DataMixingHcalWorker::addHcalSignals(const edm::Event &e) { 
00070     // fill in maps of hits
00071 
00072     LogInfo("DataMixingHcalWorker")<<"===============> adding MC signals for "<<e.id();
00073 
00074     // HBHE first
00075 
00076    Handle< HBHERecHitCollection > pHBHERecHits;
00077 
00078    const HBHERecHitCollection*  HBHERecHits = 0;
00079 
00080    if( e.getByLabel( HBHErechitCollection_.label(), pHBHERecHits) ) {
00081      HBHERecHits = pHBHERecHits.product(); // get a ptr to the product
00082      LogDebug("DataMixingHcalWorker") << "total # HBHE rechits: " << HBHERecHits->size();
00083    } 
00084    
00085  
00086    if (HBHERecHits)
00087      {
00088        // loop over rechits, storing them in a map so we can add pileup later
00089        for(HBHERecHitCollection::const_iterator it  = HBHERecHits->begin();     
00090            it != HBHERecHits->end(); ++it) {
00091 
00092          HBHERecHitStorage_.insert(HBHERecHitMap::value_type( ( it->id() ), *it ));
00093          
00094 #ifdef DEBUG     
00095          LogDebug("DataMixingHcalWorker") << "processed HBHERecHit with rawId: "
00096                                       << it->id() << "\n"
00097                                       << " rechit energy: " << it->energy();
00098 #endif
00099 
00100        }
00101      }
00102 
00103    // HO next
00104 
00105    Handle< HORecHitCollection > pHORecHits;
00106 
00107    const HORecHitCollection*  HORecHits = 0;
00108 
00109    if( e.getByLabel( HOrechitCollection_.label(), pHORecHits) ){
00110      HORecHits = pHORecHits.product(); // get a ptr to the product
00111 #ifdef DEBUG
00112      LogDebug("DataMixingHcalWorker") << "total # HO rechits: " << HORecHits->size();
00113 #endif
00114    } 
00115    
00116  
00117    if (HORecHits)
00118      {
00119        // loop over rechits, storing them in a map so we can add pileup later
00120        for(HORecHitCollection::const_iterator it  = HORecHits->begin(); 
00121            it != HORecHits->end(); ++it) {
00122 
00123          HORecHitStorage_.insert(HORecHitMap::value_type( ( it->id() ), *it ));
00124          
00125 #ifdef DEBUG     
00126          LogDebug("DataMixingHcalWorker") << "processed HORecHit with rawId: "
00127                                       << it->id() << "\n"
00128                                       << " rechit energy: " << it->energy();
00129 #endif
00130 
00131        }
00132      }
00133 
00134    // HF next
00135 
00136    Handle< HFRecHitCollection > pHFRecHits;
00137 
00138    const HFRecHitCollection*  HFRecHits = 0;
00139 
00140    if( e.getByLabel( HFrechitCollection_.label(), pHFRecHits) ) {
00141      HFRecHits = pHFRecHits.product(); // get a ptr to the product
00142 #ifdef DEBUG
00143      LogDebug("DataMixingHcalWorker") << "total # HF rechits: " << HFRecHits->size();
00144 #endif
00145    } 
00146    
00147  
00148    if (HFRecHits)
00149      {
00150        // loop over rechits, storing them in a map so we can add pileup later
00151        for(HFRecHitCollection::const_iterator it  = HFRecHits->begin(); 
00152            it != HFRecHits->end(); ++it) {
00153 
00154          HFRecHitStorage_.insert(HFRecHitMap::value_type( ( it->id() ), *it ));
00155          
00156 #ifdef DEBUG     
00157          LogDebug("DataMixingHcalWorker") << "processed HFRecHit with rawId: "
00158                                       << it->id() << "\n"
00159                                       << " rechit energy: " << it->energy();
00160 #endif
00161 
00162        }
00163      }
00164 
00165    // ZDC next
00166 
00167    Handle< ZDCRecHitCollection > pZDCRecHits;
00168 
00169    const ZDCRecHitCollection*  ZDCRecHits = 0;
00170 
00171    if( e.getByLabel( ZDCrechitCollection_.label(), pZDCRecHits) ) {
00172      ZDCRecHits = pZDCRecHits.product(); // get a ptr to the product
00173 #ifdef DEBUG
00174      LogDebug("DataMixingHcalWorker") << "total # ZDC rechits: " << ZDCRecHits->size();
00175 #endif
00176    } 
00177    
00178  
00179    if (ZDCRecHits)
00180      {
00181        // loop over rechits, storing them in a map so we can add pileup later
00182        for(ZDCRecHitCollection::const_iterator it  = ZDCRecHits->begin();       
00183            it != ZDCRecHits->end(); ++it) {
00184 
00185          ZDCRecHitStorage_.insert(ZDCRecHitMap::value_type( ( it->id() ), *it ));
00186          
00187 #ifdef DEBUG     
00188          LogDebug("DataMixingHcalWorker") << "processed ZDCRecHit with rawId: "
00189                                       << it->id() << "\n"
00190                                       << " rechit energy: " << it->energy();
00191 #endif
00192 
00193        }
00194      }
00195     
00196   } // end of addEMSignals
00197 
00198   void DataMixingHcalWorker::addHcalPileups(const int bcr, Event *e, unsigned int eventNr) {
00199   
00200     LogDebug("DataMixingHcalWorker") <<"\n===============> adding pileups from event  "<<e->id()<<" for bunchcrossing "<<bcr;
00201 
00202     // fill in maps of hits; same code as addSignals, except now applied to the pileup events
00203 
00204     // HBHE first
00205 
00206    Handle< HBHERecHitCollection > pHBHERecHits;
00207    const HBHERecHitCollection*  HBHERecHits = 0;
00208 
00209    if( e->getByLabel( HBHErechitCollection_.label(), pHBHERecHits) ) {
00210      HBHERecHits = pHBHERecHits.product(); // get a ptr to the product
00211 #ifdef DEBUG
00212      LogDebug("DataMixingHcalWorker") << "total # HEHB rechits: " << HBHERecHits->size();
00213 #endif
00214    } 
00215    
00216  
00217    if (HBHERecHits)
00218      {
00219        // loop over rechits, adding these to the existing maps
00220        for(HBHERecHitCollection::const_iterator it  = HBHERecHits->begin();
00221            it != HBHERecHits->end(); ++it) {
00222 
00223          HBHERecHitStorage_.insert(HBHERecHitMap::value_type( (it->id()), *it ));
00224          
00225 #ifdef DEBUG     
00226          LogDebug("DataMixingHcalWorker") << "processed HBHERecHit with rawId: "
00227                                       << it->id() << "\n"
00228                                       << " rechit energy: " << it->energy();
00229 #endif
00230        }
00231      }
00232     // HO Next
00233 
00234    Handle< HORecHitCollection > pHORecHits;
00235    const HORecHitCollection*  HORecHits = 0;
00236 
00237    if( e->getByLabel( HOrechitCollection_.label(), pHORecHits) ) {
00238      HORecHits = pHORecHits.product(); // get a ptr to the product
00239 #ifdef DEBUG
00240      LogDebug("DataMixingHcalWorker") << "total # HO rechits: " << HORecHits->size();
00241 #endif
00242    } 
00243    
00244  
00245    if (HORecHits)
00246      {
00247        // loop over rechits, adding these to the existing maps
00248        for(HORecHitCollection::const_iterator it  = HORecHits->begin();
00249            it != HORecHits->end(); ++it) {
00250 
00251          HORecHitStorage_.insert(HORecHitMap::value_type( (it->id()), *it ));
00252          
00253 #ifdef DEBUG     
00254          LogDebug("DataMixingHcalWorker") << "processed HORecHit with rawId: "
00255                                       << it->id() << "\n"
00256                                       << " rechit energy: " << it->energy();
00257 #endif
00258        }
00259      }
00260 
00261     // HF Next
00262 
00263    Handle< HFRecHitCollection > pHFRecHits;
00264    const HFRecHitCollection*  HFRecHits = 0;
00265 
00266    if( e->getByLabel( HFrechitCollection_.label(), pHFRecHits) ) {
00267      HFRecHits = pHFRecHits.product(); // get a ptr to the product
00268 #ifdef DEBUG
00269      LogDebug("DataMixingHcalWorker") << "total # HF rechits: " << HFRecHits->size();
00270 #endif
00271    } 
00272    
00273  
00274    if (HFRecHits)
00275      {
00276        // loop over rechits, adding these to the existing maps
00277        for(HFRecHitCollection::const_iterator it  = HFRecHits->begin();
00278            it != HFRecHits->end(); ++it) {
00279 
00280          HFRecHitStorage_.insert(HFRecHitMap::value_type( (it->id()), *it ));
00281          
00282 #ifdef DEBUG     
00283          LogDebug("DataMixingHcalWorker") << "processed HFRecHit with rawId: "
00284                                       << it->id() << "\n"
00285                                       << " rechit energy: " << it->energy();
00286 #endif
00287        }
00288      }
00289 
00290     // ZDC Next
00291 
00292    Handle< ZDCRecHitCollection > pZDCRecHits;
00293    const ZDCRecHitCollection*  ZDCRecHits = 0;
00294 
00295    if( e->getByLabel( ZDCrechitCollection_.label(), pZDCRecHits) ) {
00296      ZDCRecHits = pZDCRecHits.product(); // get a ptr to the product
00297 #ifdef DEBUG
00298      LogDebug("DataMixingHcalWorker") << "total # ZDC rechits: " << ZDCRecHits->size();
00299 #endif
00300    } 
00301    
00302  
00303    if (ZDCRecHits)
00304      {
00305        // loop over rechits, adding these to the existing maps
00306        for(ZDCRecHitCollection::const_iterator it  = ZDCRecHits->begin();
00307            it != ZDCRecHits->end(); ++it) {
00308 
00309          ZDCRecHitStorage_.insert(ZDCRecHitMap::value_type( (it->id()), *it ));
00310          
00311 #ifdef DEBUG     
00312          LogDebug("DataMixingHcalWorker") << "processed ZDCRecHit with rawId: "
00313                                       << it->id() << "\n"
00314                                       << " rechit energy: " << it->energy();
00315 #endif
00316        }
00317      }
00318 
00319   }
00320  
00321   void DataMixingHcalWorker::putHcal(edm::Event &e) {
00322 
00323     // collection of rechits to put in the event
00324     std::auto_ptr< HBHERecHitCollection > HBHErechits( new HBHERecHitCollection );
00325     std::auto_ptr< HORecHitCollection > HOrechits( new HORecHitCollection );
00326     std::auto_ptr< HFRecHitCollection > HFrechits( new HFRecHitCollection );
00327     std::auto_ptr< ZDCRecHitCollection > ZDCrechits( new ZDCRecHitCollection );
00328 
00329     // loop over the maps we have, re-making individual hits or digis if necessary.
00330     DetId formerID = 0;
00331     DetId currentID;
00332     float ESum = 0.;
00333     float HBTime = 0.;
00334     HBHERecHit OldHit;
00335     int nmatch=0;
00336 
00337     // HB first...
00338 
00339     HBHERecHitMap::const_iterator iHBchk;
00340 
00341     for(HBHERecHitMap::const_iterator iHB  = HBHERecHitStorage_.begin();
00342         iHB != HBHERecHitStorage_.end(); ++iHB) {
00343 
00344       currentID = iHB->first; 
00345 
00346       if (currentID == formerID) { // we have to add these rechits together
00347         nmatch++;                  // use this to avoid using the "count" function
00348         ESum+=(iHB->second).energy();          // on every element...
00349 
00350         iHBchk = iHB;
00351         if((iHBchk++) == HBHERecHitStorage_.end()) {  //make sure not to lose the last one
00352           HBHERecHit aHit(formerID, ESum, HBTime);
00353           HBHErechits->push_back( aHit );         
00354           // reset energy sum, nmatch
00355           ESum = 0 ;
00356           nmatch=0;       
00357         }
00358       }
00359       else {
00360         if(nmatch>0) {
00361           HBHERecHit aHit(formerID, ESum, HBTime);
00362           HBHErechits->push_back( aHit );         
00363           // reset energy sum, nmatch
00364           ESum = 0 ;
00365           nmatch=0;       
00366         }
00367         else {
00368           if(formerID>0) HBHErechits->push_back( OldHit );
00369         }
00370         
00371         iHBchk = iHB;
00372         if((iHBchk++) == HBHERecHitStorage_.end()) {  //make sure not to lose the last one
00373           HBHErechits->push_back( iHB->second );
00374         }
00375 
00376         // save pointers for next iteration
00377         OldHit = iHB->second;
00378         formerID = currentID;
00379         ESum = (iHB->second).energy();
00380         HBTime = (iHB->second).time();  // take time of first hit in sequence - is this ok?
00381       }
00382     }
00383 
00384     // HO next...
00385 
00386     // loop over the maps we have, re-making individual hits or digis if necessary.
00387     formerID = 0;
00388     ESum = 0.;
00389     float HOTime = 0.;
00390     HORecHit HOOldHit;
00391     nmatch=0;
00392 
00393     HORecHitMap::const_iterator iHOchk;
00394 
00395     for(HORecHitMap::const_iterator iHO  = HORecHitStorage_.begin();
00396         iHO != HORecHitStorage_.end(); ++iHO) {
00397 
00398       currentID = iHO->first; 
00399 
00400       if (currentID == formerID) { // we have to add these rechits together
00401         nmatch++;                  // use this to avoid using the "count" function
00402         ESum+=(iHO->second).energy();          // on every element...
00403 
00404         iHOchk = iHO;
00405         if((iHOchk++) == HORecHitStorage_.end()) {  //make sure not to lose the last one
00406           HORecHit aHit(formerID, ESum, HOTime);
00407           HOrechits->push_back( aHit );   
00408           // reset energy sum, nmatch
00409           ESum = 0 ;
00410           nmatch=0;       
00411         }
00412       }
00413       else {
00414         if(nmatch>0) {
00415           HORecHit aHit(formerID, ESum, HOTime);
00416           HOrechits->push_back( aHit );   
00417           // reset energy sum, nmatch
00418           ESum = 0 ;
00419           nmatch=0;       
00420         }
00421         else {
00422           if(formerID>0) HOrechits->push_back( HOOldHit );
00423         }
00424         
00425         iHOchk = iHO;
00426         if((iHOchk++) == HORecHitStorage_.end()) {  //make sure not to lose the last one
00427           HOrechits->push_back( iHO->second );
00428         }
00429 
00430         // save pointers for next iteration
00431         HOOldHit = iHO->second;
00432         formerID = currentID;
00433         ESum = (iHO->second).energy();
00434         HOTime = (iHO->second).time();  // take time of first hit in sequence - is this ok?
00435       }
00436     }
00437 
00438     // HF next...
00439 
00440     // loop over the maps we have, re-making individual hits or digis if necessary.
00441     formerID = 0;
00442     ESum = 0.;
00443     float HFTime = 0.;
00444     HFRecHit HFOldHit;
00445     nmatch=0;
00446 
00447     HFRecHitMap::const_iterator iHFchk;
00448 
00449     for(HFRecHitMap::const_iterator iHF  = HFRecHitStorage_.begin();
00450         iHF != HFRecHitStorage_.end(); ++iHF) {
00451 
00452       currentID = iHF->first; 
00453 
00454       if (currentID == formerID) { // we have to add these rechits together
00455         nmatch++;                  // use this to avoid using the "count" function
00456         ESum+=(iHF->second).energy();          // on every element...
00457 
00458         iHFchk = iHF;
00459         if((iHFchk++) == HFRecHitStorage_.end()) {  //make sure not to lose the last one
00460           HFRecHit aHit(formerID, ESum, HFTime);
00461           HFrechits->push_back( aHit );   
00462           // reset energy sum, nmatch
00463           ESum = 0 ;
00464           nmatch=0;       
00465         }
00466       }
00467       else {
00468         if(nmatch>0) {
00469           HFRecHit aHit(formerID, ESum, HFTime);
00470           HFrechits->push_back( aHit );   
00471           // reset energy sum, nmatch
00472           ESum = 0 ;
00473           nmatch=0;       
00474         }
00475         else {
00476           if(formerID>0) HFrechits->push_back( HFOldHit );
00477         }
00478         
00479         iHFchk = iHF;
00480         if((iHFchk++) == HFRecHitStorage_.end()) {  //make sure not to lose the last one
00481           HFrechits->push_back( iHF->second );
00482         }
00483 
00484         // save pointers for next iteration
00485         HFOldHit = iHF->second;
00486         formerID = currentID;
00487         ESum = (iHF->second).energy();
00488         HFTime = (iHF->second).time();  // take time of first hit in sequence - is this ok?
00489       }
00490     }
00491 
00492     // ZDC next...
00493 
00494     // loop over the maps we have, re-making individual hits or digis if necessary.
00495     formerID = 0;
00496     ESum = 0.;
00497     float ZDCTime = 0.;
00498     ZDCRecHit ZOldHit;
00499     nmatch=0;
00500 
00501     ZDCRecHitMap::const_iterator iZDCchk;
00502 
00503     for(ZDCRecHitMap::const_iterator iZDC  = ZDCRecHitStorage_.begin();
00504         iZDC != ZDCRecHitStorage_.end(); ++iZDC) {
00505 
00506       currentID = iZDC->first; 
00507 
00508       if (currentID == formerID) { // we have to add these rechits together
00509         nmatch++;                  // use this to avoid using the "count" function
00510         ESum+=(iZDC->second).energy();          // on every element...
00511 
00512         iZDCchk = iZDC;
00513         if((iZDCchk++) == ZDCRecHitStorage_.end()) {  //make sure not to lose the last one
00514           ZDCRecHit aHit(formerID, ESum, ZDCTime);
00515           ZDCrechits->push_back( aHit );          
00516           // reset energy sum, nmatch
00517           ESum = 0 ;
00518           nmatch=0;       
00519         }
00520       }
00521       else {
00522         if(nmatch>0) {
00523           ZDCRecHit aHit(formerID, ESum, ZDCTime);
00524           ZDCrechits->push_back( aHit );          
00525           // reset energy sum, nmatch
00526           ESum = 0 ;
00527           nmatch=0;       
00528         }
00529         else {
00530           if(formerID>0) ZDCrechits->push_back( ZOldHit );
00531         }
00532         
00533         iZDCchk = iZDC;
00534         if((iZDCchk++) == ZDCRecHitStorage_.end()) {  //make sure not to lose the last one
00535           ZDCrechits->push_back( iZDC->second );
00536         }
00537 
00538         // save pointers for next iteration
00539         ZOldHit = iZDC->second;
00540         formerID = currentID;
00541         ESum = (iZDC->second).energy();
00542         ZDCTime = (iZDC->second).time();  // take time of first hit in sequence - is this ok?
00543       }
00544     }
00545 
00546     // put the collection of recunstructed hits in the event   
00547     LogInfo("DataMixingHcalWorker") << "total # HBHE Merged rechits: " << HBHErechits->size() ;
00548     LogInfo("DataMixingHcalWorker") << "total # HO Merged rechits: " << HOrechits->size() ;
00549     LogInfo("DataMixingHcalWorker") << "total # HF Merged rechits: " << HFrechits->size() ;
00550     LogInfo("DataMixingHcalWorker") << "total # ZDC Merged rechits: " << ZDCrechits->size() ;
00551 
00552     e.put( HBHErechits, HBHERecHitCollectionDM_ );
00553     e.put( HOrechits, HORecHitCollectionDM_ );
00554     e.put( HFrechits, HFRecHitCollectionDM_ );
00555     e.put( ZDCrechits, ZDCRecHitCollectionDM_ );
00556 
00557     // clear local storage after this event
00558     HBHERecHitStorage_.clear();
00559     HORecHitStorage_.clear();
00560     HFRecHitStorage_.clear();
00561     ZDCRecHitStorage_.clear();
00562 
00563   }
00564 
00565 } //edm

Generated on Tue Jun 9 17:47:18 2009 for CMSSW by  doxygen 1.5.4