CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HiEvtPlaneFlatCalib.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HiEvtPlaneFlatCalib
4 // Class: HiEvtPlaneFlatCalib
5 //
14 //
15 // Original Author: Stephen Sanders
16 // Created: Sat Jun 26 16:04:04 EDT 2010
17 // $Id: HiEvtPlaneFlatCalib.cc,v 1.7 2012/02/15 11:04:09 eulisse Exp $
18 //
19 //
20 
21 
22 // system include files
23 #include <memory>
24 
25 // user include files
28 
31 
34 
35 #include "Math/Vector3D.h"
36 
40 
43 
49 
56 
58 #include "TFile.h"
59 #include "TH1.h"
60 #include "TH2D.h"
61 #include "TH2F.h"
62 #include "TTree.h"
63 #include "TH1I.h"
64 #include "TF1.h"
65 #include "TList.h"
66 #include "TString.h"
67 #include <time.h>
68 #include <cstdlib>
69 #include <vector>
70 
72 
73 using namespace std;
74 using namespace hi;
75 
76 //
77 // class declaration
78 //
79 
81  public:
82  explicit HiEvtPlaneFlatCalib(const edm::ParameterSet&);
84 
85  private:
86  virtual void beginJob() ;
87  virtual void analyze(const edm::Event&, const edm::EventSetup&);
88  virtual void endJob() ;
89 
90  // ----------member data ---------------------------
92  // const CentralityBins * cbins_;
94  int vs_sell; // vertex collection size
95  float vzr_sell;
96  float vzErr_sell;
97 
98  static const int NumCentBins=9;
99  double wcent[10];
100 
101  TH1D * hcent;
102  TH1D * hvtx;
103  TH1D * flatXhist[NumEPNames];
104  TH1D * flatYhist[NumEPNames];
105  TH1D * flatCnthist[NumEPNames];
106 
107  TH1D * flatXDBhist[NumEPNames];
108  TH1D * flatYDBhist[NumEPNames];
109  TH1D * hPsi[NumEPNames];
110  TH1D * hPsiFlat[NumEPNames];
111  TH1D * hPsiFlatCent[NumEPNames][NumCentBins];
112  TH1D * hPsiFlatSub1[NumEPNames];
113  TH1D * hPsiFlatSub2[NumEPNames];
114  Double_t epang[NumEPNames];
117  int nRP;
119 
120 };
121 
122 //
123 // constants, enums and typedefs
124 //
125 
126 //
127 // static data member definitions
128 //
129 
130 //
131 // constructors and destructor
132 //
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 }
186 
187 
189 {
190 
191  // do anything here that needs to be done at desctruction time
192  // (e.g. close files, deallocate resources etc.)
193 
194 }
195 
196 
197 //
198 // member functions
199 //
200 
201 // ------------ method called to produce the data ------------
202 void
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 }
284 
285 // ------------ method called once each job just before starting event loop ------------
286 void
288 {
289 }
290 
291 // ------------ method called once each job just after ending the event loop ------------
292 void
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 }
302 
303 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
const int EPOrder[]
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const std::string EPNames[]
void beginJob()
Definition: Breakpoints.cc:15
int iEvent
Definition: GenABIO.cc:243
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual void analyze(const edm::Event &, const edm::EventSetup &)
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
const T & get() const
Definition: EventSetup.h:55
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
T * make() const
make new ROOT object
int RPNameIndx[50]
Definition: RPFlatParams.h:14
static const int NumEPNames
HiEvtPlaneFlatCalib(const edm::ParameterSet &)
CentralityProvider * centrality_
edm::Service< TFileService > fs