CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes

HiEvtPlaneFlatCalib Class Reference

#include <HiEvtPlaneFlatten/HiEvtPlaneFlatCalib/src/HiEvtPlaneFlatCalib.cc>

Inheritance diagram for HiEvtPlaneFlatCalib:
edm::EDAnalyzer

List of all members.

Public Member Functions

 HiEvtPlaneFlatCalib (const edm::ParameterSet &)
 ~HiEvtPlaneFlatCalib ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()

Private Attributes

CentralityProvidercentrality_
Double_t epang [NumEPNames]
HiEvtPlaneFlattenflat [NumEPNames]
TH1D * flatCnthist [NumEPNames]
TH1D * flatXDBhist [NumEPNames]
TH1D * flatXhist [NumEPNames]
TH1D * flatYDBhist [NumEPNames]
TH1D * flatYhist [NumEPNames]
edm::Service< TFileServicefs
bool genFlatPsi_
TH1D * hcent
TH1D * hPsi [NumEPNames]
TH1D * hPsiFlat [NumEPNames]
TH1D * hPsiFlatCent [NumEPNames][NumCentBins]
TH1D * hPsiFlatSub1 [NumEPNames]
TH1D * hPsiFlatSub2 [NumEPNames]
TH1D * hvtx
int nRP
RPFlatParamsrpFlat
int vs_sell
float vzErr_sell
float vzr_sell
double wcent [10]

Static Private Attributes

static const int NumCentBins = 9

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 79 of file HiEvtPlaneFlatCalib.cc.


Constructor & Destructor Documentation

HiEvtPlaneFlatCalib::HiEvtPlaneFlatCalib ( const edm::ParameterSet iConfig) [explicit]

Definition at line 132 of file HiEvtPlaneFlatCalib.cc.

References data, hi::EPNames, hi::EPOrder, edm::ParameterSet::getUntrackedParameter(), i, j, TFileDirectory::make(), TFileDirectory::mkdir(), pileupCalc::nbins, and hi::NumEPNames.

{
  genFlatPsi_ = iConfig.getUntrackedParameter<bool>("genFlatPsi_",true);

  //  NumCentBins=9;
  wcent[0] = 0;
  wcent[1] = 5;
  wcent[2] = 10;
  wcent[3] = 20;
  wcent[4] = 30;
  wcent[5] = 40;
  wcent[6] = 50;
  wcent[7] = 60;
  wcent[8] = 70;
  wcent[9] = 100;

  //now do what ever other initialization is needed
  //  cbins_ = 0;
  centrality_ = 0;
  hcent = fs->make<TH1D>("cent","cent",41,0,40);
  hvtx = fs->make<TH1D>("vtx","vtx",1000,-50,50);
  //setting before 8:35 CDT on 11Jan2011
  Int_t FlatOrder = 21;
  for(int i = 0; i<NumEPNames; i++) {
    TFileDirectory subdir = fs->mkdir(Form("%s",EPNames[i].data()));
    flat[i] = new HiEvtPlaneFlatten();
    flat[i]->Init(FlatOrder,11,4,EPNames[i],EPOrder[i]);
    int nbins = flat[i]->GetHBins();
    flatXhist[i] = subdir.make<TH1D>(Form("x_%s",EPNames[i].data()),Form("x_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
    flatYhist[i] = subdir.make<TH1D>(Form("y_%s",EPNames[i].data()),Form("y_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
    flatCnthist[i] = subdir.make<TH1D>(Form("cnt_%s",EPNames[i].data()),Form("cnt_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
    Double_t psirange = 4;
    if(EPOrder[i]==2 ) psirange = 2;
    if(EPOrder[i]==3 ) psirange = 1.5;
    if(EPOrder[i]==4 ) psirange = 1;
    if(EPOrder[i]==5) psirange = 0.8;
    if(EPOrder[i]==6) psirange = 0.6;
    hPsi[i] = subdir.make<TH1D>("psi","psi",800,-psirange,psirange);
    hPsi[i]->SetXTitle("#Psi");
    hPsi[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
    hPsiFlat[i] = subdir.make<TH1D>("psiFlat","psiFlat",800,-psirange,psirange);
    hPsiFlat[i]->SetXTitle("#Psi");
    hPsiFlat[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
    for(int j = 0; j<NumCentBins; j++) {
      TString hname = Form("psiFlat_%d_%d",(int) wcent[j],(int) wcent[j+1]);
      hPsiFlatCent[i][j] = subdir.make<TH1D>(hname.Data(),hname.Data(),800,-psirange,psirange);
      hPsiFlatCent[i][j]->SetXTitle("#Psi");
      hPsiFlatCent[i][j]->SetYTitle(Form("Counts (%d<cent#leq%d%c)",(int) wcent[j],(int) wcent[j+1],'%'));
    }  

  }
  
}
HiEvtPlaneFlatCalib::~HiEvtPlaneFlatCalib ( )

Definition at line 187 of file HiEvtPlaneFlatCalib.cc.

{
 
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

void HiEvtPlaneFlatCalib::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 202 of file HiEvtPlaneFlatCalib.cc.

References newFWLiteAna::bin, compare_using_db::compare, hi::EPNames, edm::EventSetup::get(), edm::Event::getByLabel(), i, j, hi::NumEPNames, dt_dqm_sourceclient_common_cff::reco, RPFlatParams::EP::RPNameIndx, RPFlatParams::EP::x, and RPFlatParams::EP::y.

{
  using namespace edm;
  using namespace reco;
  //
  //Get Centrality
  //
  if(!centrality_) centrality_ = new CentralityProvider(iSetup);

   centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
   int bin = centrality_->getBin();
  double centval = 2.5*bin+1.25;
  hcent->Fill(bin);
  //
  //Get Vertex
  //
  edm::Handle<reco::VertexCollection> vertexCollection3;
  iEvent.getByLabel("hiSelectedVertex",vertexCollection3);
  const reco::VertexCollection * vertices3 = vertexCollection3.product();
  vs_sell = vertices3->size();
  if(vs_sell>0) {
    vzr_sell = vertices3->begin()->z();
    vzErr_sell = vertices3->begin()->zError();
  } else
    vzr_sell = -999.9;
  //
  //Get Flattening Parameters
  //
  if(genFlatPsi_) {
    edm::ESHandle<RPFlatParams> flatparmsDB_;
    iSetup.get<HeavyIonRPRcd>().get(flatparmsDB_);
    int flatTableSize = flatparmsDB_->m_table.size();
    for(int i = 0; i<flatTableSize; i++) {
      const RPFlatParams::EP* thisBin = &(flatparmsDB_->m_table[i]);
      for(int j = 0; j<NumEPNames; j++) {
        int indx = thisBin->RPNameIndx[j];
        if(indx>=0) {
            flat[indx]->SetXDB(i, thisBin->x[j]);
            flat[indx]->SetYDB(i, thisBin->y[j]);
        }
      }
    }
  }
    
  //
  //Get Event Planes
  //

  Handle<reco::EvtPlaneCollection> evtPlanes;
  iEvent.getByLabel("hiEvtPlane","recoLevel",evtPlanes);

  if(!evtPlanes.isValid()){
    //cout << "Error! Can't get hiEvtPlane product!" << endl;
    return ;
  }

  for (EvtPlaneCollection::const_iterator rp = evtPlanes->begin();rp !=evtPlanes->end(); rp++) {
    if(rp->angle() > -5) {
      string baseName = rp->label();
      for(int i = 0; i< NumEPNames; i++) {
        if(EPNames[i].compare(baseName)==0) {
          double psiFlat = flat[i]->GetFlatPsi(rp->angle(),vzr_sell,bin);
          epang[i]=psiFlat;
          if(EPNames[i].compare(rp->label())==0) {
            flat[i]->Fill(rp->angle(),vzr_sell,bin);
            if(i==0)  hvtx->Fill(vzr_sell);

            if(centval<=80) hPsi[i]->Fill(rp->angle());
            if(genFlatPsi_) {
              if(centval<=80) hPsiFlat[i]->Fill(psiFlat);
              for(int j = 0; j<NumCentBins; j++) {
                if(centval>wcent[j]&&centval<=wcent[j+1]) hPsiFlatCent[i][j]->Fill(psiFlat);
              }
            }
          } 
        }
      }
    }    
  }

}
void HiEvtPlaneFlatCalib::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 286 of file HiEvtPlaneFlatCalib.cc.

{
}
void HiEvtPlaneFlatCalib::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 292 of file HiEvtPlaneFlatCalib.cc.

References i, j, and hi::NumEPNames.

                            {
  for(int i = 0; i<NumEPNames; i++) {
    for(int j = 0; j<flat[i]->GetHBins();j++) {
      flatXhist[i]->SetBinContent(j+1,flat[i]->GetX(j));
      flatYhist[i]->SetBinContent(j+1,flat[i]->GetY(j));
      flatCnthist[i]->SetBinContent(j+1,flat[i]->GetCnt(j));
    }
  }
}

Member Data Documentation

Definition at line 92 of file HiEvtPlaneFlatCalib.cc.

Double_t HiEvtPlaneFlatCalib::epang[NumEPNames] [private]

Definition at line 113 of file HiEvtPlaneFlatCalib.cc.

Definition at line 114 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatCnthist[NumEPNames] [private]

Definition at line 104 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatXDBhist[NumEPNames] [private]

Definition at line 106 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatXhist[NumEPNames] [private]

Definition at line 102 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatYDBhist[NumEPNames] [private]

Definition at line 107 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatYhist[NumEPNames] [private]

Definition at line 103 of file HiEvtPlaneFlatCalib.cc.

Definition at line 90 of file HiEvtPlaneFlatCalib.cc.

Definition at line 117 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hcent [private]

Definition at line 100 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsi[NumEPNames] [private]

Definition at line 108 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlat[NumEPNames] [private]

Definition at line 109 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlatCent[NumEPNames][NumCentBins] [private]

Definition at line 110 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlatSub1[NumEPNames] [private]

Definition at line 111 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlatSub2[NumEPNames] [private]

Definition at line 112 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hvtx [private]

Definition at line 101 of file HiEvtPlaneFlatCalib.cc.

int HiEvtPlaneFlatCalib::nRP [private]

Definition at line 116 of file HiEvtPlaneFlatCalib.cc.

const int HiEvtPlaneFlatCalib::NumCentBins = 9 [static, private]

Definition at line 97 of file HiEvtPlaneFlatCalib.cc.

Definition at line 115 of file HiEvtPlaneFlatCalib.cc.

Definition at line 93 of file HiEvtPlaneFlatCalib.cc.

Definition at line 95 of file HiEvtPlaneFlatCalib.cc.

Definition at line 94 of file HiEvtPlaneFlatCalib.cc.

double HiEvtPlaneFlatCalib::wcent[10] [private]

Definition at line 98 of file HiEvtPlaneFlatCalib.cc.