Go to the documentation of this file.00001 #include "Calibration/HcalAlCaRecoProducers/interface/AlCaDiJetsProducer.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00004 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00005 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00006 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00007 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00008 #include "FWCore/Utilities/interface/Exception.h"
00009
00010 using namespace edm;
00011 using namespace std;
00012 using namespace reco;
00013
00014 namespace cms
00015 {
00016
00017 AlCaDiJetsProducer::AlCaDiJetsProducer(const edm::ParameterSet& iConfig)
00018 {
00019 jetsInput_ = iConfig.getParameter<edm::InputTag>("jetsInput");
00020 ecalLabels_=iConfig.getParameter<std::vector<edm::InputTag> >("ecalInputs");
00021 hbheInput_ = iConfig.getParameter<edm::InputTag>("hbheInput");
00022 hoInput_ = iConfig.getParameter<edm::InputTag>("hoInput");
00023 hfInput_ = iConfig.getParameter<edm::InputTag>("hfInput");
00024 allowMissingInputs_ = true;
00025
00026
00027 produces<CaloJetCollection>("DiJetsBackToBackCollection");
00028 produces<EcalRecHitCollection>("DiJetsEcalRecHitCollection");
00029 produces<HBHERecHitCollection>("DiJetsHBHERecHitCollection");
00030 produces<HORecHitCollection>("DiJetsHORecHitCollection");
00031 produces<HFRecHitCollection>("DiJetsHFRecHitCollection");
00032
00033 }
00034 void AlCaDiJetsProducer::beginJob()
00035 {
00036 }
00037
00038 AlCaDiJetsProducer::~AlCaDiJetsProducer()
00039 {
00040
00041
00042 }
00043
00044
00045
00046 void
00047 AlCaDiJetsProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00048 {
00049
00050 double pi = 4.*atan(1.);
00051
00052 std::auto_ptr<CaloJetCollection> result (new CaloJetCollection);
00053 std::auto_ptr<EcalRecHitCollection> miniDiJetsEcalRecHitCollection(new EcalRecHitCollection);
00054 std::auto_ptr<HBHERecHitCollection> miniDiJetsHBHERecHitCollection(new HBHERecHitCollection);
00055 std::auto_ptr<HORecHitCollection> miniDiJetsHORecHitCollection(new HORecHitCollection);
00056 std::auto_ptr<HFRecHitCollection> miniDiJetsHFRecHitCollection(new HFRecHitCollection);
00057
00058 edm::ESHandle<CaloGeometry> pG;
00059 iSetup.get<CaloGeometryRecord>().get(pG);
00060 const CaloGeometry* geo = pG.product();
00061
00062
00063
00064 vector<CaloJet> jetv;
00065
00066 CaloJet fJet1, fJet2, fJet3;
00067 edm::Handle<CaloJetCollection> jets;
00068 iEvent.getByLabel(jetsInput_, jets);
00069 int iflag_select = 0;
00070 if(jets->size()>1){
00071 fJet1 = (*jets)[0];
00072 fJet2 = (*jets)[1];
00073 double dphi = fabs(fJet1.phi() - fJet2.phi());
00074 if(dphi > pi){dphi = 2*pi - dphi;}
00075 double degreedphi = dphi*180./pi;
00076 if(fabs(degreedphi-180)<30.){iflag_select = 1;}
00077 }
00078 if(iflag_select == 1){
00079 result->push_back(fJet1);
00080 result->push_back(fJet2);
00081 jetv.push_back(fJet1);
00082 jetv.push_back(fJet2);
00083 if(jets->size()>2){
00084 fJet3 = (*jets)[2];
00085 result->push_back(fJet3);
00086 jetv.push_back(fJet3);
00087 }
00088 } else {
00089 iEvent.put( result, "DiJetsBackToBackCollection");
00090 iEvent.put( miniDiJetsEcalRecHitCollection,"DiJetsEcalRecHitCollection");
00091 iEvent.put( miniDiJetsHBHERecHitCollection, "DiJetsHBHERecHitCollection");
00092 iEvent.put( miniDiJetsHORecHitCollection, "DiJetsHORecHitCollection");
00093 iEvent.put( miniDiJetsHFRecHitCollection, "DiJetsHFRecHitCollection");
00094 return;
00095 }
00096
00097
00098
00099 std::vector<edm::InputTag>::const_iterator i;
00100 for (i=ecalLabels_.begin(); i!=ecalLabels_.end(); i++) {
00101 edm::Handle<EcalRecHitCollection> ec;
00102 iEvent.getByLabel(*i,ec);
00103 for(EcalRecHitCollection::const_iterator ecItr = (*ec).begin();
00104 ecItr != (*ec).end(); ++ecItr)
00105 {
00106
00107 GlobalPoint pos = geo->getPosition(ecItr->detid());
00108 double phihit = pos.phi();
00109 double etahit = pos.eta();
00110 int iflag_select = 0;
00111 for(unsigned int i=0; i<jetv.size(); i++){
00112 double deta = fabs(etahit - jetv[i].eta());
00113 double dphi = fabs(phihit - jetv[i].phi());
00114 if(dphi > pi) dphi = 2*pi - dphi;
00115 double dr = sqrt(deta*deta+dphi*dphi);
00116 if(dr < 1.4){iflag_select = 1;}
00117 }
00118 if(iflag_select==1){miniDiJetsEcalRecHitCollection->push_back(*ecItr);}
00119 }
00120
00121 }
00122
00123
00124
00125 edm::Handle<HBHERecHitCollection> hbhe;
00126 iEvent.getByLabel(hbheInput_,hbhe);
00127 for(HBHERecHitCollection::const_iterator hbheItr=hbhe->begin();
00128 hbheItr!=hbhe->end(); hbheItr++)
00129 {
00130 GlobalPoint pos = geo->getPosition(hbheItr->detid());
00131 double phihit = pos.phi();
00132 double etahit = pos.eta();
00133 int iflag_select = 0;
00134 for(unsigned int i=0; i<jetv.size(); i++){
00135 double deta = fabs(etahit - jetv[i].eta());
00136 double dphi = fabs(phihit - jetv[i].phi());
00137 if(dphi > pi) dphi = 2*pi - dphi;
00138 double dr = sqrt(deta*deta+dphi*dphi);
00139 if(dr < 1.4){iflag_select = 1;}
00140 }
00141 if(iflag_select==1){miniDiJetsHBHERecHitCollection->push_back(*hbheItr);}
00142 }
00143
00144
00145
00146
00147
00148 edm::Handle<HORecHitCollection> ho;
00149 iEvent.getByLabel(hoInput_,ho);
00150 for(HORecHitCollection::const_iterator hoItr=ho->begin();
00151 hoItr!=ho->end(); hoItr++)
00152 {
00153 GlobalPoint pos = geo->getPosition(hoItr->detid());
00154 double phihit = pos.phi();
00155 double etahit = pos.eta();
00156 int iflag_select = 0;
00157 for(unsigned int i=0; i<jetv.size(); i++){
00158 double deta = fabs(etahit - jetv[i].eta());
00159 double dphi = fabs(phihit - jetv[i].phi());
00160 if(dphi > pi) dphi = 2*pi - dphi;
00161 double dr = sqrt(deta*deta+dphi*dphi);
00162 if(dr < 1.4){iflag_select = 1;}
00163 }
00164 if(iflag_select==1){miniDiJetsHORecHitCollection->push_back(*hoItr);}
00165 }
00166
00167
00168
00169
00170 edm::Handle<HFRecHitCollection> hf;
00171 iEvent.getByLabel(hfInput_,hf);
00172 for(HFRecHitCollection::const_iterator hfItr=hf->begin();
00173 hfItr!=hf->end(); hfItr++)
00174 {
00175 GlobalPoint pos = geo->getPosition(hfItr->detid());
00176 double phihit = pos.phi();
00177 double etahit = pos.eta();
00178 int iflag_select = 0;
00179 for(unsigned int i=0; i<jetv.size(); i++){
00180 double deta = fabs(etahit - jetv[i].eta());
00181 double dphi = fabs(phihit - jetv[i].phi());
00182 if(dphi > pi) dphi = 2*pi - dphi;
00183 double dr = sqrt(deta*deta+dphi*dphi);
00184 if(dr < 1.4){iflag_select = 1;}
00185 }
00186 if(iflag_select==1){miniDiJetsHFRecHitCollection->push_back(*hfItr);}
00187 }
00188
00189
00190
00191
00192 iEvent.put( result, "DiJetsBackToBackCollection");
00193 iEvent.put( miniDiJetsEcalRecHitCollection,"DiJetsEcalRecHitCollection");
00194 iEvent.put( miniDiJetsHBHERecHitCollection, "DiJetsHBHERecHitCollection");
00195 iEvent.put( miniDiJetsHORecHitCollection, "DiJetsHORecHitCollection");
00196 iEvent.put( miniDiJetsHFRecHitCollection, "DiJetsHFRecHitCollection");
00197 }
00198 }