Go to the documentation of this file.00001
00002 #include "JetMETCorrections/MinBias/interface/MinBias.h"
00003 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00004 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00005
00006 using namespace std;
00007 namespace cms
00008 {
00009 MinBias::MinBias(const edm::ParameterSet& iConfig)
00010 {
00011
00012 hbheLabel_= iConfig.getParameter<std::string>("hbheInput");
00013 hoLabel_=iConfig.getParameter<std::string>("hoInput");
00014 hfLabel_=iConfig.getParameter<std::string>("hfInput");
00015 allowMissingInputs_=iConfig.getUntrackedParameter<bool>("AllowMissingInputs",false);
00016
00017 fOutputFileName = iConfig.getUntrackedParameter<std::string>("HistOutFile");
00018
00019 }
00020
00021
00022 MinBias::~MinBias()
00023 {
00024
00025
00026
00027
00028 }
00029
00030 void MinBias::beginJob()
00031 {
00032 hOutputFile = new TFile( fOutputFileName.c_str(), "RECREATE" ) ;
00033 myTree = new TTree("RecJet","RecJet Tree");
00034 myTree->Branch("mydet", &mydet, "mydet/I");
00035 myTree->Branch("mysubd", &mysubd, "mysubd/I");
00036 myTree->Branch("depth", &depth, "depth/I");
00037 myTree->Branch("ieta", &ieta, "ieta/I");
00038 myTree->Branch("iphi", &iphi, "iphi/I");
00039 myTree->Branch("eta", &eta, "eta/F");
00040 myTree->Branch("phi", &phi, "phi/F");
00041 myTree->Branch("mom1", &mom1, "mom1/F");
00042 myTree->Branch("mom2", &mom2, "mom2/F");
00043 myTree->Branch("mom3", &mom3, "mom3/F");
00044 myTree->Branch("mom4", &mom4, "mom4/F");
00045
00046 ievent = 0;
00047
00048 }
00049
00050 void MinBias::endJob()
00051 {
00052
00053 const std::vector<DetId>& did = geo->getSubdetectorGeometry( DetId::Hcal, 1 )->getValidDetIds() ;
00054 int i=0;
00055 for(std::vector<DetId>::const_iterator id=did.begin(); id != did.end(); id++)
00056 {
00057
00058 GlobalPoint pos = geo->getPosition(*id);
00059 mydet = ((*id).rawId()>>28)&0xF;
00060 mysubd = ((*id).rawId()>>25)&0x7;
00061 depth = HcalDetId(*id).depth();
00062 ieta = HcalDetId(*id).ieta();
00063 iphi = HcalDetId(*id).iphi();
00064 phi = pos.phi();
00065 eta = pos.eta();
00066 if ( theFillDetMap0[*id] > 0. )
00067 {
00068 mom1 = theFillDetMap1[*id]/theFillDetMap0[*id];
00069 mom2 = theFillDetMap2[*id]/theFillDetMap0[*id]-(mom1*mom1);
00070 mom3 = theFillDetMap3[*id]/theFillDetMap0[*id]-3.*mom1*theFillDetMap2[*id]/theFillDetMap0[*id]+
00071 2.*pow(mom2,3);
00072 mom4 = (theFillDetMap4[*id]-4.*mom1*theFillDetMap3[*id]+6.*pow(mom1,2)*theFillDetMap2[*id])/theFillDetMap0[*id]-3.*pow(mom1,4);
00073
00074 } else
00075 {
00076 mom1 = 0.; mom2 = 0.; mom3 = 0.; mom4 = 0.;
00077 }
00078 std::cout<<" Detector "<<(*id).rawId()<<" mydet "<<mydet<<" "<<mysubd<<" "<<depth<<" "<<
00079 HcalDetId(*id).subdet()<<" "<<ieta<<" "<<iphi<<" "<<pos.eta()<<" "<<pos.phi()<<std::endl;
00080 std::cout<<" Energy "<<mom1<<" "<<mom2<<std::endl;
00081 myTree->Fill();
00082 i++;
00083
00084 }
00085 std::cout<<" The number of CaloDet records "<<did.size()<<std::endl;
00086 std::cout<<" The number of Hcal records "<<i<<std::endl;
00087
00088
00089 std::cout << "===== Start writing user histograms =====" << std::endl;
00090 hOutputFile->SetCompressionLevel(2);
00091 hOutputFile->cd();
00092 myTree->Write();
00093 hOutputFile->Close() ;
00094 std::cout << "===== End writing user histograms =======" << std::endl;
00095 }
00096
00097
00098
00099
00100
00101
00102
00103 void
00104 MinBias::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00105 {
00106
00107 using namespace edm;
00108 if(ievent == 0 ){
00109 edm::ESHandle<CaloGeometry> pG;
00110 iSetup.get<CaloGeometryRecord>().get(pG);
00111 geo = pG.product();
00112 std::vector<DetId> did = geo->getValidDetIds();
00113
00114 for(std::vector<DetId>::const_iterator id=did.begin(); id != did.end(); id++)
00115 {
00116 if( (*id).det() == DetId::Hcal ) {
00117 theFillDetMap0[*id] = 0.;
00118 theFillDetMap1[*id] = 0.;
00119 theFillDetMap2[*id] = 0.;
00120 theFillDetMap3[*id] = 0.;
00121 theFillDetMap4[*id] = 0.;
00122 }
00123 }
00124 }
00125
00126
00127
00128 if (!hbheLabel_.empty()) {
00129 edm::Handle<HBHERecHitCollection> hbhe;
00130 iEvent.getByLabel(hbheLabel_,hbhe);
00131 if (!hbhe.isValid()) {
00132
00133 if (!allowMissingInputs_) {
00134 *hbhe;
00135 }
00136 } else {
00137 for(HBHERecHitCollection::const_iterator hbheItr = (*hbhe).begin();
00138 hbheItr != (*hbhe).end(); ++hbheItr)
00139 {
00140 DetId id = (hbheItr)->detid();
00141 if( (*hbheItr).energy() > 0. ) std::cout<<" Energy = "<<(*hbheItr).energy()<<std::endl;
00142 theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
00143 theFillDetMap1[id] = theFillDetMap1[id]+(*hbheItr).energy();
00144 theFillDetMap2[id] = theFillDetMap2[id]+pow((*hbheItr).energy(),2);
00145 theFillDetMap3[id] = theFillDetMap3[id]+pow((*hbheItr).energy(),3);
00146 theFillDetMap4[id] = theFillDetMap4[id]+pow((*hbheItr).energy(),4);
00147 }
00148 }
00149 }
00150
00151 if (!hoLabel_.empty()) {
00152 edm::Handle<HORecHitCollection> ho;
00153 iEvent.getByLabel(hoLabel_,ho);
00154 if (!ho.isValid()) {
00155
00156 if (!allowMissingInputs_) {
00157 *ho;
00158 }
00159 } else {
00160 for(HORecHitCollection::const_iterator hoItr = (*ho).begin();
00161 hoItr != (*ho).end(); ++hoItr)
00162 {
00163 DetId id = (hoItr)->detid();
00164 theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
00165 theFillDetMap1[id] = theFillDetMap1[id]+(*hoItr).energy();
00166 theFillDetMap2[id] = theFillDetMap2[id]+pow((*hoItr).energy(),2);
00167 theFillDetMap3[id] = theFillDetMap3[id]+pow((*hoItr).energy(),3);
00168 theFillDetMap4[id] = theFillDetMap4[id]+pow((*hoItr).energy(),4);
00169 }
00170 }
00171 }
00172
00173 if (!hfLabel_.empty()) {
00174 edm::Handle<HFRecHitCollection> hf;
00175 iEvent.getByLabel(hfLabel_,hf);
00176 if (!hf.isValid()) {
00177
00178 if (!allowMissingInputs_) {
00179 *hf;
00180 }
00181 } else {
00182 for(HFRecHitCollection::const_iterator hfItr = (*hf).begin();
00183 hfItr != (*hf).end(); ++hfItr)
00184 {
00185 DetId id = (hfItr)->detid();
00186 theFillDetMap0[id] = theFillDetMap0[id]+ 1.;
00187 theFillDetMap1[id] = theFillDetMap1[id]+(*hfItr).energy();
00188 theFillDetMap2[id] = theFillDetMap2[id]+pow((*hfItr).energy(),2);
00189 theFillDetMap3[id] = theFillDetMap3[id]+pow((*hfItr).energy(),3);
00190 theFillDetMap4[id] = theFillDetMap4[id]+pow((*hfItr).energy(),4);
00191 }
00192 }
00193 }
00194
00195 }
00196 }