CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoHI/HiEvtPlaneAlgos/src/HiEvtPlaneFlatCalib.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HiEvtPlaneFlatCalib
00004 // Class:      HiEvtPlaneFlatCalib
00005 // 
00014 //
00015 // Original Author:  Stephen Sanders
00016 //         Created:  Sat Jun 26 16:04:04 EDT 2010
00017 // $Id: HiEvtPlaneFlatCalib.cc,v 1.7 2012/02/15 11:04:09 eulisse Exp $
00018 //
00019 //
00020 
00021 
00022 // system include files
00023 #include <memory>
00024 
00025 // user include files
00026 #include "FWCore/Framework/interface/Frameworkfwd.h"
00027 #include "FWCore/Framework/interface/EDAnalyzer.h"
00028 
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 
00032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00033 #include "DataFormats/HeavyIonEvent/interface/CentralityProvider.h"
00034 #include "Math/Vector3D.h"
00035 
00036 #include "DataFormats/Common/interface/Handle.h"
00037 #include "FWCore/Framework/interface/ESHandle.h"
00038 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00039 
00040 #include "DataFormats/HeavyIonEvent/interface/Centrality.h"
00041 #include "DataFormats/HeavyIonEvent/interface/CentralityBins.h"
00042 
00043 #include "DataFormats/HeavyIonEvent/interface/EvtPlane.h"
00044 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00045 #include "DataFormats/TrackReco/interface/Track.h"
00046 #include "DataFormats/VertexReco/interface/Vertex.h"
00047 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00048 
00049 #include "FWCore/ServiceRegistry/interface/Service.h"
00050 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00051 #include "CondFormats/DataRecord/interface/HeavyIonRPRcd.h"
00052 #include "CondFormats/HIObjects/interface/CentralityTable.h"
00053 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00054 #include "CondFormats/HIObjects/interface/RPFlatParams.h"
00055 
00056 #include "RecoHI/HiEvtPlaneAlgos/interface/HiEvtPlaneFlatten.h"
00057 #include "TFile.h"
00058 #include "TH1.h"
00059 #include "TH2D.h"
00060 #include "TH2F.h"
00061 #include "TTree.h"
00062 #include "TH1I.h"
00063 #include "TF1.h"
00064 #include "TList.h"
00065 #include "TString.h"
00066 #include <time.h>
00067 #include <cstdlib>
00068 #include <vector>
00069 
00070 #include "RecoHI/HiEvtPlaneAlgos/interface/HiEvtPlaneList.h"
00071 
00072 using namespace std;
00073 using namespace hi;
00074 
00075 //
00076 // class declaration
00077 //
00078 
00079 class HiEvtPlaneFlatCalib : public edm::EDAnalyzer {
00080    public:
00081       explicit HiEvtPlaneFlatCalib(const edm::ParameterSet&);
00082       ~HiEvtPlaneFlatCalib();
00083 
00084    private:
00085       virtual void beginJob() ;
00086       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00087       virtual void endJob() ;
00088       
00089       // ----------member data ---------------------------
00090   edm::Service<TFileService> fs;
00091   //  const CentralityBins * cbins_;
00092   CentralityProvider * centrality_;
00093   int vs_sell;   // vertex collection size
00094   float vzr_sell;
00095   float vzErr_sell;
00096 
00097   static const  int NumCentBins=9;
00098   double wcent[10];
00099 
00100   TH1D * hcent;
00101   TH1D * hvtx;
00102   TH1D * flatXhist[NumEPNames];
00103   TH1D * flatYhist[NumEPNames];
00104   TH1D * flatCnthist[NumEPNames];
00105 
00106   TH1D * flatXDBhist[NumEPNames];
00107   TH1D * flatYDBhist[NumEPNames];
00108   TH1D * hPsi[NumEPNames];
00109   TH1D * hPsiFlat[NumEPNames];
00110   TH1D * hPsiFlatCent[NumEPNames][NumCentBins];
00111   TH1D * hPsiFlatSub1[NumEPNames];
00112   TH1D * hPsiFlatSub2[NumEPNames];
00113   Double_t epang[NumEPNames];
00114   HiEvtPlaneFlatten * flat[NumEPNames];
00115   RPFlatParams * rpFlat;
00116   int nRP;
00117   bool genFlatPsi_;
00118 
00119 };
00120 
00121 //
00122 // constants, enums and typedefs
00123 //
00124 
00125 //
00126 // static data member definitions
00127 //
00128 
00129 //
00130 // constructors and destructor
00131 //
00132 HiEvtPlaneFlatCalib::HiEvtPlaneFlatCalib(const edm::ParameterSet& iConfig)
00133 {
00134   genFlatPsi_ = iConfig.getUntrackedParameter<bool>("genFlatPsi_",true);
00135 
00136   //  NumCentBins=9;
00137   wcent[0] = 0;
00138   wcent[1] = 5;
00139   wcent[2] = 10;
00140   wcent[3] = 20;
00141   wcent[4] = 30;
00142   wcent[5] = 40;
00143   wcent[6] = 50;
00144   wcent[7] = 60;
00145   wcent[8] = 70;
00146   wcent[9] = 100;
00147 
00148   //now do what ever other initialization is needed
00149   //  cbins_ = 0;
00150   centrality_ = 0;
00151   hcent = fs->make<TH1D>("cent","cent",41,0,40);
00152   hvtx = fs->make<TH1D>("vtx","vtx",1000,-50,50);
00153   //setting before 8:35 CDT on 11Jan2011
00154   Int_t FlatOrder = 21;
00155   for(int i = 0; i<NumEPNames; i++) {
00156     TFileDirectory subdir = fs->mkdir(Form("%s",EPNames[i].data()));
00157     flat[i] = new HiEvtPlaneFlatten();
00158     flat[i]->Init(FlatOrder,11,4,EPNames[i],EPOrder[i]);
00159     int nbins = flat[i]->GetHBins();
00160     flatXhist[i] = subdir.make<TH1D>(Form("x_%s",EPNames[i].data()),Form("x_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
00161     flatYhist[i] = subdir.make<TH1D>(Form("y_%s",EPNames[i].data()),Form("y_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
00162     flatCnthist[i] = subdir.make<TH1D>(Form("cnt_%s",EPNames[i].data()),Form("cnt_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
00163     Double_t psirange = 4;
00164     if(EPOrder[i]==2 ) psirange = 2;
00165     if(EPOrder[i]==3 ) psirange = 1.5;
00166     if(EPOrder[i]==4 ) psirange = 1;
00167     if(EPOrder[i]==5) psirange = 0.8;
00168     if(EPOrder[i]==6) psirange = 0.6;
00169     hPsi[i] = subdir.make<TH1D>("psi","psi",800,-psirange,psirange);
00170     hPsi[i]->SetXTitle("#Psi");
00171     hPsi[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
00172     hPsiFlat[i] = subdir.make<TH1D>("psiFlat","psiFlat",800,-psirange,psirange);
00173     hPsiFlat[i]->SetXTitle("#Psi");
00174     hPsiFlat[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
00175     for(int j = 0; j<NumCentBins; j++) {
00176       TString hname = Form("psiFlat_%d_%d",(int) wcent[j],(int) wcent[j+1]);
00177       hPsiFlatCent[i][j] = subdir.make<TH1D>(hname.Data(),hname.Data(),800,-psirange,psirange);
00178       hPsiFlatCent[i][j]->SetXTitle("#Psi");
00179       hPsiFlatCent[i][j]->SetYTitle(Form("Counts (%d<cent#leq%d%c)",(int) wcent[j],(int) wcent[j+1],'%'));
00180     }  
00181 
00182   }
00183   
00184 }
00185 
00186 
00187 HiEvtPlaneFlatCalib::~HiEvtPlaneFlatCalib()
00188 {
00189  
00190    // do anything here that needs to be done at desctruction time
00191    // (e.g. close files, deallocate resources etc.)
00192 
00193 }
00194 
00195 
00196 //
00197 // member functions
00198 //
00199 
00200 // ------------ method called to produce the data  ------------
00201 void
00202 HiEvtPlaneFlatCalib::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00203 {
00204   using namespace edm;
00205   using namespace reco;
00206   //
00207   //Get Centrality
00208   //
00209   if(!centrality_) centrality_ = new CentralityProvider(iSetup);
00210 
00211    centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
00212    int bin = centrality_->getBin();
00213   double centval = 2.5*bin+1.25;
00214   hcent->Fill(bin);
00215   //
00216   //Get Vertex
00217   //
00218   edm::Handle<reco::VertexCollection> vertexCollection3;
00219   iEvent.getByLabel("hiSelectedVertex",vertexCollection3);
00220   const reco::VertexCollection * vertices3 = vertexCollection3.product();
00221   vs_sell = vertices3->size();
00222   if(vs_sell>0) {
00223     vzr_sell = vertices3->begin()->z();
00224     vzErr_sell = vertices3->begin()->zError();
00225   } else
00226     vzr_sell = -999.9;
00227   //
00228   //Get Flattening Parameters
00229   //
00230   if(genFlatPsi_) {
00231     edm::ESHandle<RPFlatParams> flatparmsDB_;
00232     iSetup.get<HeavyIonRPRcd>().get(flatparmsDB_);
00233     int flatTableSize = flatparmsDB_->m_table.size();
00234     for(int i = 0; i<flatTableSize; i++) {
00235       const RPFlatParams::EP* thisBin = &(flatparmsDB_->m_table[i]);
00236       for(int j = 0; j<NumEPNames; j++) {
00237         int indx = thisBin->RPNameIndx[j];
00238         if(indx>=0) {
00239             flat[indx]->SetXDB(i, thisBin->x[j]);
00240             flat[indx]->SetYDB(i, thisBin->y[j]);
00241         }
00242       }
00243     }
00244   }
00245     
00246   //
00247   //Get Event Planes
00248   //
00249 
00250   Handle<reco::EvtPlaneCollection> evtPlanes;
00251   iEvent.getByLabel("hiEvtPlane","recoLevel",evtPlanes);
00252 
00253   if(!evtPlanes.isValid()){
00254     //cout << "Error! Can't get hiEvtPlane product!" << endl;
00255     return ;
00256   }
00257 
00258   for (EvtPlaneCollection::const_iterator rp = evtPlanes->begin();rp !=evtPlanes->end(); rp++) {
00259     if(rp->angle() > -5) {
00260       string baseName = rp->label();
00261       for(int i = 0; i< NumEPNames; i++) {
00262         if(EPNames[i].compare(baseName)==0) {
00263           double psiFlat = flat[i]->GetFlatPsi(rp->angle(),vzr_sell,bin);
00264           epang[i]=psiFlat;
00265           if(EPNames[i].compare(rp->label())==0) {
00266             flat[i]->Fill(rp->angle(),vzr_sell,bin);
00267             if(i==0)  hvtx->Fill(vzr_sell);
00268 
00269             if(centval<=80) hPsi[i]->Fill(rp->angle());
00270             if(genFlatPsi_) {
00271               if(centval<=80) hPsiFlat[i]->Fill(psiFlat);
00272               for(int j = 0; j<NumCentBins; j++) {
00273                 if(centval>wcent[j]&&centval<=wcent[j+1]) hPsiFlatCent[i][j]->Fill(psiFlat);
00274               }
00275             }
00276           } 
00277         }
00278       }
00279     }    
00280   }
00281 
00282 }
00283 
00284 // ------------ method called once each job just before starting event loop  ------------
00285 void 
00286 HiEvtPlaneFlatCalib::beginJob()
00287 {
00288 }
00289 
00290 // ------------ method called once each job just after ending the event loop  ------------
00291 void 
00292 HiEvtPlaneFlatCalib::endJob() {
00293   for(int i = 0; i<NumEPNames; i++) {
00294     for(int j = 0; j<flat[i]->GetHBins();j++) {
00295       flatXhist[i]->SetBinContent(j+1,flat[i]->GetX(j));
00296       flatYhist[i]->SetBinContent(j+1,flat[i]->GetY(j));
00297       flatCnthist[i]->SetBinContent(j+1,flat[i]->GetCnt(j));
00298     }
00299   }
00300 }
00301 
00302 //define this as a plug-in
00303 DEFINE_FWK_MODULE(HiEvtPlaneFlatCalib);