CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/Calibration/HcalAlCaRecoProducers/src/AlCaHcalNoiseProducer.cc

Go to the documentation of this file.
00001 /*
00002 Original author Grigory Safronov
00003 
00004 27/03/09 - compilation from :
00005 HLTrigger/special/src/HLTHcalNoiseFilter.cc
00006 Calibration/HcalAlCaRecoProducers/src/AlCaEcalHcalReadoutsProducer.cc
00007 Calibration/HcalIsolatedTrackReco/src/SubdetFEDSelector.cc
00008 
00009 */
00010 
00011 
00012 #include "Calibration/HcalAlCaRecoProducers/interface/AlCaHcalNoiseProducer.h"
00013 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00014 
00015 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00016 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00017 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00018 
00019 #include "EventFilter/RawDataCollector/interface/RawDataFEDSelector.h"
00020 
00021 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00022 #include "DataFormats/CaloTowers/interface/CaloTower.h"
00023 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00024 #include "DataFormats/JetReco/interface/CaloJet.h"
00025 #include "DataFormats/METReco/interface/CaloMET.h"
00026 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00027 #include "DataFormats/Math/interface/deltaR.h"
00028 
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 
00031 
00032 AlCaHcalNoiseProducer::AlCaHcalNoiseProducer(const edm::ParameterSet& iConfig)
00033 {
00034    JetSource_ = iConfig.getParameter<edm::InputTag>("JetSource");
00035    MetSource_ = iConfig.getParameter<edm::InputTag>("MetSource");
00036    TowerSource_ = iConfig.getParameter<edm::InputTag>("TowerSource");
00037    useMet_ = iConfig.getParameter<bool>("UseMET");
00038    useJet_ = iConfig.getParameter<bool>("UseJet");
00039    MetCut_ = iConfig.getParameter<double>("MetCut");
00040    JetMinE_ = iConfig.getParameter<double>("JetMinE");
00041    JetHCALminEnergyFraction_ = iConfig.getParameter<double>("JetHCALminEnergyFraction");
00042 
00043    hoLabel_ = iConfig.getParameter<edm::InputTag>("hoInput");
00044    hfLabel_ = iConfig.getParameter<edm::InputTag>("hfInput");
00045    hbheLabel_ = iConfig.getParameter<edm::InputTag>("hbheInput");
00046    ecalLabels_= iConfig.getParameter<std::vector<edm::InputTag> >("ecalInputs");
00047    ecalPSLabel_=iConfig.getParameter<edm::InputTag>("ecalPSInput");
00048    rawInLabel_=iConfig.getParameter<edm::InputTag>("rawInput");
00049 
00050   
00051    //register products
00052    produces<HBHERecHitCollection>("HBHERecHitCollectionFHN");
00053    produces<HORecHitCollection>("HORecHitCollectionFHN");
00054    produces<HFRecHitCollection>("HFRecHitCollectionFHN");
00055    
00056    produces<EcalRecHitCollection>("EcalRecHitCollectionFHN");
00057    produces<EcalRecHitCollection>("PSEcalRecHitCollectionFHN");
00058     
00059    produces<FEDRawDataCollection>("HcalFEDsFHN");
00060 
00061 }
00062 
00063 
00064 AlCaHcalNoiseProducer::~AlCaHcalNoiseProducer()
00065 {
00066 }
00067 
00068 
00069 // ------------ method called to produce the data  ------------
00070 void
00071 AlCaHcalNoiseProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00072 {
00073    bool acceptEvent=false;
00074 
00075   // filtering basing on HLTrigger/special/src/HLTHcalNoiseFilter.cc:
00076 
00077    bool isAnomalous_BasedOnMET = false;
00078    bool isAnomalous_BasedOnEnergyFraction=false; 
00079 
00080    if (useMet_)
00081      {
00082        edm::Handle <reco::CaloMETCollection> metHandle;
00083        iEvent.getByLabel(MetSource_, metHandle);
00084        const reco::CaloMETCollection *metCol = metHandle.product();
00085        const reco::CaloMET met = metCol->front();
00086     
00087        if(met.pt() > MetCut_) isAnomalous_BasedOnMET=true;
00088      }
00089        
00090    if (useJet_)
00091      {
00092        edm::Handle<reco::CaloJetCollection> calojetHandle;
00093        iEvent.getByLabel(JetSource_,calojetHandle);
00094        
00095        edm::Handle<CaloTowerCollection> towerHandle;
00096        iEvent.getByLabel(TowerSource_, towerHandle);
00097 
00098        std::vector<CaloTower> TowerContainer;
00099        std::vector<reco::CaloJet> JetContainer;
00100        TowerContainer.clear();
00101        JetContainer.clear();
00102        CaloTower seedTower;
00103        nEvents++;
00104        for(reco::CaloJetCollection::const_iterator calojetIter = calojetHandle->begin();calojetIter != calojetHandle->end();++calojetIter) {
00105          if( ((calojetIter->et())*cosh(calojetIter->eta()) > JetMinE_) && (calojetIter->energyFractionHadronic() > JetHCALminEnergyFraction_) ) {
00106            JetContainer.push_back(*calojetIter);
00107            double maxTowerE = 0.0;
00108            for(CaloTowerCollection::const_iterator kal = towerHandle->begin(); kal != towerHandle->end(); kal++) {
00109              double dR = deltaR((*calojetIter).eta(),(*calojetIter).phi(),(*kal).eta(),(*kal).phi());
00110              if( (dR < 0.50) && (kal->p() > maxTowerE) ) {
00111                maxTowerE = kal->p();
00112                seedTower = *kal;
00113              }
00114            }
00115            TowerContainer.push_back(seedTower);
00116          }
00117          
00118        }
00119        if(JetContainer.size() > 0) {
00120          nAnomalousEvents++;
00121          isAnomalous_BasedOnEnergyFraction = true;
00122        }
00123      }
00124    
00125    acceptEvent=((useMet_&&isAnomalous_BasedOnMET)||(useJet_&&isAnomalous_BasedOnEnergyFraction));
00126 
00128 
00129     //Create empty output collections
00130    
00131   std::auto_ptr<HBHERecHitCollection> miniHBHERecHitCollection(new HBHERecHitCollection);
00132   std::auto_ptr<HORecHitCollection> miniHORecHitCollection(new HORecHitCollection);
00133   std::auto_ptr<HFRecHitCollection> miniHFRecHitCollection(new HFRecHitCollection);
00134 
00135   std::auto_ptr<EcalRecHitCollection> outputEColl(new EcalRecHitCollection);
00136   std::auto_ptr<EcalRecHitCollection> outputESColl(new EcalRecHitCollection);
00137 
00138   std::auto_ptr<FEDRawDataCollection> outputFEDs(new FEDRawDataCollection);
00139 
00140 
00141   // if good event get and save all colletions
00142   if (acceptEvent)
00143     {
00144       edm::Handle<HBHERecHitCollection> hbhe;
00145       edm::Handle<HORecHitCollection> ho;
00146       edm::Handle<HFRecHitCollection> hf;
00147       
00148       iEvent.getByLabel(hbheLabel_,hbhe);
00149       iEvent.getByLabel(hoLabel_,ho);
00150       iEvent.getByLabel(hfLabel_,hf);
00151 
00152       edm::Handle<EcalRecHitCollection> pRecHits;
00153       iEvent.getByLabel(ecalPSLabel_,pRecHits);
00154 
00155       // temporary collection of EB+EE recHits
00156 
00157       std::auto_ptr<EcalRecHitCollection> tmpEcalRecHitCollection(new EcalRecHitCollection);
00158       
00159       std::vector<edm::InputTag>::const_iterator i;
00160       for (i=ecalLabels_.begin(); i!=ecalLabels_.end(); i++) 
00161         {
00162           edm::Handle<EcalRecHitCollection> ec;
00163           iEvent.getByLabel(*i,ec);
00164           for(EcalRecHitCollection::const_iterator recHit = (*ec).begin(); recHit != (*ec).end(); ++recHit)
00165             {
00166               tmpEcalRecHitCollection->push_back(*recHit);
00167             }
00168         }    
00169 
00171 
00173       const HBHERecHitCollection Hithbhe = *(hbhe.product());
00174       for(HBHERecHitCollection::const_iterator hbheItr=Hithbhe.begin(); hbheItr!=Hithbhe.end(); hbheItr++)
00175         {
00176           miniHBHERecHitCollection->push_back(*hbheItr);
00177         }
00178       const HORecHitCollection Hitho = *(ho.product());
00179       for(HORecHitCollection::const_iterator hoItr=Hitho.begin(); hoItr!=Hitho.end(); hoItr++)
00180         {
00181           miniHORecHitCollection->push_back(*hoItr);
00182         }
00183       
00184       const HFRecHitCollection Hithf = *(hf.product());
00185       for(HFRecHitCollection::const_iterator hfItr=Hithf.begin(); hfItr!=Hithf.end(); hfItr++)
00186         {
00187           miniHFRecHitCollection->push_back(*hfItr);
00188         }
00190 
00192       for (std::vector<EcalRecHit>::const_iterator ehit=tmpEcalRecHitCollection->begin(); ehit!=tmpEcalRecHitCollection->end(); ehit++) 
00193         {
00194           outputEColl->push_back(*ehit);
00195         }
00197 
00198       // write PS
00199       const EcalRecHitCollection& psrechits = *(pRecHits.product());
00200       
00201       for(EcalRecHitCollection::const_iterator i=psrechits.begin(); i!=psrechits.end(); i++) 
00202         {
00203           outputESColl->push_back( *i );
00204         }
00205       
00206 
00207       // get HCAL FEDs
00208       edm::Handle<FEDRawDataCollection> rawIn;
00209       iEvent.getByLabel(rawInLabel_,rawIn);
00210       
00211       std::vector<int> selFEDs;
00212       for (int i=FEDNumbering::MINHCALFEDID; i<=FEDNumbering::MAXHCALFEDID; i++)
00213         {
00214           selFEDs.push_back(i);
00215         }
00217 
00218       // Copying FEDs :
00219       const FEDRawDataCollection *rdc=rawIn.product();
00220       
00221       //   if ( ( rawData[i].provenance()->processName() != e.processHistory().rbegin()->processName() ) )
00222       //       continue ; // skip all raw collections not produced by the current process
00223       
00224       for ( int j=0; j< FEDNumbering::MAXFEDID; ++j ) 
00225         {
00226           bool rightFED=false;
00227           for (uint32_t k=0; k<selFEDs.size(); k++)
00228             {
00229               if (j==selFEDs[k])
00230                 {
00231                   rightFED=true;
00232                 }
00233             }
00234           if (!rightFED) continue;
00235           const FEDRawData & fedData = rdc->FEDData(j);
00236           size_t size=fedData.size();
00237           
00238           if ( size > 0 ) 
00239             {
00240               // this fed has data -- lets copy it
00241               FEDRawData & fedDataProd = outputFEDs->FEDData(j);
00242               if ( fedDataProd.size() != 0 ) {
00243                 //         std::cout << " More than one FEDRawDataCollection with data in FED ";
00244                 //         std::cout << j << " Skipping the 2nd\n";
00245                 continue;
00246               }
00247               fedDataProd.resize(size);
00248               unsigned char *dataProd=fedDataProd.data();
00249               const unsigned char *data=fedData.data();
00250               for ( unsigned int k=0; k<size; ++k ) {
00251                 dataProd[k]=data[k];
00252               }
00253             }
00254         }
00256     }
00257 
00258   //Put selected information in the event
00259   iEvent.put( miniHBHERecHitCollection, "HBHERecHitCollectionFHN");
00260   iEvent.put( miniHORecHitCollection, "HORecHitCollectionFHN");
00261   iEvent.put( miniHFRecHitCollection, "HFRecHitCollectionFHN");
00262   iEvent.put( outputEColl, "EcalRecHitCollectionFHN"); 
00263   iEvent.put( outputESColl, "PSEcalRecHitCollectionFHN");
00264   iEvent.put( outputFEDs, "HcalFEDsFHN"); 
00265 }