CMS 3D CMS Logo

MinBias.cc

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

Generated on Tue Jun 9 17:39:36 2009 for CMSSW by  doxygen 1.5.4