CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | 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

Public Member Functions

 HiEvtPlaneFlatCalib (const edm::ParameterSet &)
 
 ~HiEvtPlaneFlatCalib ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

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
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 80 of file HiEvtPlaneFlatCalib.cc.

Constructor & Destructor Documentation

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

Definition at line 133 of file HiEvtPlaneFlatCalib.cc.

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

134 {
135  genFlatPsi_ = iConfig.getUntrackedParameter<bool>("genFlatPsi_",true);
136 
137  // NumCentBins=9;
138  wcent[0] = 0;
139  wcent[1] = 5;
140  wcent[2] = 10;
141  wcent[3] = 20;
142  wcent[4] = 30;
143  wcent[5] = 40;
144  wcent[6] = 50;
145  wcent[7] = 60;
146  wcent[8] = 70;
147  wcent[9] = 100;
148 
149  //now do what ever other initialization is needed
150  // cbins_ = 0;
151  centrality_ = 0;
152  hcent = fs->make<TH1D>("cent","cent",41,0,40);
153  hvtx = fs->make<TH1D>("vtx","vtx",1000,-50,50);
154  //setting before 8:35 CDT on 11Jan2011
155  Int_t FlatOrder = 21;
156  for(int i = 0; i<NumEPNames; i++) {
157  TFileDirectory subdir = fs->mkdir(Form("%s",EPNames[i].data()));
158  flat[i] = new HiEvtPlaneFlatten();
159  flat[i]->Init(FlatOrder,11,4,EPNames[i],EPOrder[i]);
160  int nbins = flat[i]->GetHBins();
161  flatXhist[i] = subdir.make<TH1D>(Form("x_%s",EPNames[i].data()),Form("x_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
162  flatYhist[i] = subdir.make<TH1D>(Form("y_%s",EPNames[i].data()),Form("y_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
163  flatCnthist[i] = subdir.make<TH1D>(Form("cnt_%s",EPNames[i].data()),Form("cnt_%s",EPNames[i].data()),nbins,-0.5,nbins-0.5);
164  Double_t psirange = 4;
165  if(EPOrder[i]==2 ) psirange = 2;
166  if(EPOrder[i]==3 ) psirange = 1.5;
167  if(EPOrder[i]==4 ) psirange = 1;
168  if(EPOrder[i]==5) psirange = 0.8;
169  if(EPOrder[i]==6) psirange = 0.6;
170  hPsi[i] = subdir.make<TH1D>("psi","psi",800,-psirange,psirange);
171  hPsi[i]->SetXTitle("#Psi");
172  hPsi[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
173  hPsiFlat[i] = subdir.make<TH1D>("psiFlat","psiFlat",800,-psirange,psirange);
174  hPsiFlat[i]->SetXTitle("#Psi");
175  hPsiFlat[i]->SetYTitle(Form("Counts (cent<80%c)",'%'));
176  for(int j = 0; j<NumCentBins; j++) {
177  TString hname = Form("psiFlat_%d_%d",(int) wcent[j],(int) wcent[j+1]);
178  hPsiFlatCent[i][j] = subdir.make<TH1D>(hname.Data(),hname.Data(),800,-psirange,psirange);
179  hPsiFlatCent[i][j]->SetXTitle("#Psi");
180  hPsiFlatCent[i][j]->SetYTitle(Form("Counts (%d<cent#leq%d%c)",(int) wcent[j],(int) wcent[j+1],'%'));
181  }
182 
183  }
184 
185 }
T getUntrackedParameter(std::string const &, T const &) const
const int EPOrder[]
int i
Definition: DBlmapReader.cc:9
static const int NumCentBins
TH1D * flatYhist[NumEPNames]
void Init(int order, int ncentbins, const int centbinCompression, std::string tag, int vord)
const std::string EPNames[]
HiEvtPlaneFlatten * flat[NumEPNames]
TH1D * flatXhist[NumEPNames]
TH1D * flatCnthist[NumEPNames]
int j
Definition: DBlmapReader.cc:9
TH1D * hPsi[NumEPNames]
TH1D * hPsiFlat[NumEPNames]
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
T * make() const
make new ROOT object
static const int NumEPNames
TH1D * hPsiFlatCent[NumEPNames][NumCentBins]
CentralityProvider * centrality_
edm::Service< TFileService > fs
HiEvtPlaneFlatCalib::~HiEvtPlaneFlatCalib ( )

Definition at line 188 of file HiEvtPlaneFlatCalib.cc.

189 {
190 
191  // do anything here that needs to be done at desctruction time
192  // (e.g. close files, deallocate resources etc.)
193 
194 }

Member Function Documentation

void HiEvtPlaneFlatCalib::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 203 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, hitfit::return, RPFlatParams::EP::RPNameIndx, RPFlatParams::EP::x, and RPFlatParams::EP::y.

204 {
205  using namespace edm;
206  using namespace reco;
207  //
208  //Get Centrality
209  //
210  if(!centrality_) centrality_ = new CentralityProvider(iSetup);
211 
212  centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
213  int bin = centrality_->getBin();
214  double centval = 2.5*bin+1.25;
215  hcent->Fill(bin);
216  //
217  //Get Vertex
218  //
219  edm::Handle<reco::VertexCollection> vertexCollection3;
220  iEvent.getByLabel("hiSelectedVertex",vertexCollection3);
221  const reco::VertexCollection * vertices3 = vertexCollection3.product();
222  vs_sell = vertices3->size();
223  if(vs_sell>0) {
224  vzr_sell = vertices3->begin()->z();
225  vzErr_sell = vertices3->begin()->zError();
226  } else
227  vzr_sell = -999.9;
228  //
229  //Get Flattening Parameters
230  //
231  if(genFlatPsi_) {
232  edm::ESHandle<RPFlatParams> flatparmsDB_;
233  iSetup.get<HeavyIonRPRcd>().get(flatparmsDB_);
234  int flatTableSize = flatparmsDB_->m_table.size();
235  for(int i = 0; i<flatTableSize; i++) {
236  const RPFlatParams::EP* thisBin = &(flatparmsDB_->m_table[i]);
237  for(int j = 0; j<NumEPNames; j++) {
238  int indx = thisBin->RPNameIndx[j];
239  if(indx>=0) {
240  flat[indx]->SetXDB(i, thisBin->x[j]);
241  flat[indx]->SetYDB(i, thisBin->y[j]);
242  }
243  }
244  }
245  }
246 
247  //
248  //Get Event Planes
249  //
250 
252  iEvent.getByLabel("hiEvtPlane","recoLevel",evtPlanes);
253 
254  if(!evtPlanes.isValid()){
255  //cout << "Error! Can't get hiEvtPlane product!" << endl;
256  return ;
257  }
258 
259  for (EvtPlaneCollection::const_iterator rp = evtPlanes->begin();rp !=evtPlanes->end(); rp++) {
260  if(rp->angle() > -5) {
261  string baseName = rp->label();
262  for(int i = 0; i< NumEPNames; i++) {
263  if(EPNames[i].compare(baseName)==0) {
264  double psiFlat = flat[i]->GetFlatPsi(rp->angle(),vzr_sell,bin);
265  epang[i]=psiFlat;
266  if(EPNames[i].compare(rp->label())==0) {
267  flat[i]->Fill(rp->angle(),vzr_sell,bin);
268  if(i==0) hvtx->Fill(vzr_sell);
269 
270  if(centval<=80) hPsi[i]->Fill(rp->angle());
271  if(genFlatPsi_) {
272  if(centval<=80) hPsiFlat[i]->Fill(psiFlat);
273  for(int j = 0; j<NumCentBins; j++) {
274  if(centval>wcent[j]&&centval<=wcent[j+1]) hPsiFlatCent[i][j]->Fill(psiFlat);
275  }
276  }
277  }
278  }
279  }
280  }
281  }
282 
283 }
int i
Definition: DBlmapReader.cc:9
static const int NumCentBins
double GetFlatPsi(double psi, double vtx, double cent)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const std::string EPNames[]
HiEvtPlaneFlatten * flat[NumEPNames]
void SetXDB(int indx, double val)
void SetYDB(int indx, double val)
void newEvent(const edm::Event &ev, const edm::EventSetup &iSetup)
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
TH1D * hPsi[NumEPNames]
TH1D * hPsiFlat[NumEPNames]
void Fill(double psi, double vtx, int centbin)
const T & get() const
Definition: EventSetup.h:55
int RPNameIndx[50]
Definition: RPFlatParams.h:14
static const int NumEPNames
Double_t epang[NumEPNames]
TH1D * hPsiFlatCent[NumEPNames][NumCentBins]
CentralityProvider * centrality_
void HiEvtPlaneFlatCalib::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 287 of file HiEvtPlaneFlatCalib.cc.

288 {
289 }
void HiEvtPlaneFlatCalib::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 293 of file HiEvtPlaneFlatCalib.cc.

References i, j, and hi::NumEPNames.

293  {
294  for(int i = 0; i<NumEPNames; i++) {
295  for(int j = 0; j<flat[i]->GetHBins();j++) {
296  flatXhist[i]->SetBinContent(j+1,flat[i]->GetX(j));
297  flatYhist[i]->SetBinContent(j+1,flat[i]->GetY(j));
298  flatCnthist[i]->SetBinContent(j+1,flat[i]->GetCnt(j));
299  }
300  }
301 }
int i
Definition: DBlmapReader.cc:9
TH1D * flatYhist[NumEPNames]
HiEvtPlaneFlatten * flat[NumEPNames]
TH1D * flatXhist[NumEPNames]
TH1D * flatCnthist[NumEPNames]
int j
Definition: DBlmapReader.cc:9
static const int NumEPNames

Member Data Documentation

CentralityProvider* HiEvtPlaneFlatCalib::centrality_
private

Definition at line 93 of file HiEvtPlaneFlatCalib.cc.

Double_t HiEvtPlaneFlatCalib::epang[NumEPNames]
private

Definition at line 114 of file HiEvtPlaneFlatCalib.cc.

HiEvtPlaneFlatten* HiEvtPlaneFlatCalib::flat[NumEPNames]
private

Definition at line 115 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatCnthist[NumEPNames]
private

Definition at line 105 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatXDBhist[NumEPNames]
private

Definition at line 107 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatXhist[NumEPNames]
private

Definition at line 103 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatYDBhist[NumEPNames]
private

Definition at line 108 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::flatYhist[NumEPNames]
private

Definition at line 104 of file HiEvtPlaneFlatCalib.cc.

edm::Service<TFileService> HiEvtPlaneFlatCalib::fs
private

Definition at line 91 of file HiEvtPlaneFlatCalib.cc.

bool HiEvtPlaneFlatCalib::genFlatPsi_
private

Definition at line 118 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hcent
private

Definition at line 101 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsi[NumEPNames]
private

Definition at line 109 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlat[NumEPNames]
private

Definition at line 110 of file HiEvtPlaneFlatCalib.cc.

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

Definition at line 111 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlatSub1[NumEPNames]
private

Definition at line 112 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hPsiFlatSub2[NumEPNames]
private

Definition at line 113 of file HiEvtPlaneFlatCalib.cc.

TH1D* HiEvtPlaneFlatCalib::hvtx
private

Definition at line 102 of file HiEvtPlaneFlatCalib.cc.

int HiEvtPlaneFlatCalib::nRP
private

Definition at line 117 of file HiEvtPlaneFlatCalib.cc.

const int HiEvtPlaneFlatCalib::NumCentBins =9
staticprivate

Definition at line 98 of file HiEvtPlaneFlatCalib.cc.

RPFlatParams* HiEvtPlaneFlatCalib::rpFlat
private

Definition at line 116 of file HiEvtPlaneFlatCalib.cc.

int HiEvtPlaneFlatCalib::vs_sell
private

Definition at line 94 of file HiEvtPlaneFlatCalib.cc.

float HiEvtPlaneFlatCalib::vzErr_sell
private

Definition at line 96 of file HiEvtPlaneFlatCalib.cc.

float HiEvtPlaneFlatCalib::vzr_sell
private

Definition at line 95 of file HiEvtPlaneFlatCalib.cc.

double HiEvtPlaneFlatCalib::wcent[10]
private

Definition at line 99 of file HiEvtPlaneFlatCalib.cc.