StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFlowPhiWgtMaker.cxx
1 //
3 // $Id: StFlowPhiWgtMaker.cxx,v 1.7 2004/12/09 23:47:09 posk Exp $
4 //
5 // Authors: Art Poskanzer and Jamie Dunlop, May 2003
6 //
8 //
9 // Description: Maker to produce PhiWgt files for odd and even harmonics
10 // using StFlowEvent
11 //
13 
14 #include <Stiostream.h>
15 #include <stdlib.h>
16 #include <math.h>
17 #include "StMaker.h"
18 #include "StFlowPhiWgtMaker.h"
19 #include "StFlowMaker/StFlowMaker.h"
20 #include "StFlowMaker/StFlowEvent.h"
21 #include "StFlowMaker/StFlowConstants.h"
22 #include "StFlowMaker/StFlowSelection.h"
23 #include "PhysicalConstants.h"
24 #include "TFile.h"
25 #include "TString.h"
26 #include "TH1.h"
27 #include "TOrdCollection.h"
28 #include "StMessMgr.h"
29 #include "TText.h"
30 #define PR(x) cout << "##### FlowPhiWgt: " << (#x) << " = " << (x) << endl;
31 
32 ClassImp(StFlowPhiWgtMaker)
33 
34 //-----------------------------------------------------------------------
35 
36 StFlowPhiWgtMaker::StFlowPhiWgtMaker(const Char_t* name): StMaker(name),
37  MakerName(name) {
38  pFlowSelect = new StFlowSelection();
39 }
40 
41 //-----------------------------------------------------------------------
42 
43 StFlowPhiWgtMaker::~StFlowPhiWgtMaker() {
44 }
45 
46 //-----------------------------------------------------------------------
47 
49  // Make histograms
50 
51  // Get a pointer to StFlowEvent
52  StFlowMaker* pFlowMaker = NULL;
53  pFlowMaker = (StFlowMaker*)GetMaker("Flow");
54  if (pFlowMaker) pFlowEvent = pFlowMaker->FlowEventPointer();
55  if (pFlowEvent && pFlowSelect->Select(pFlowEvent)) { // event selected
56  FillParticleHistograms(); // fill particle histograms
57  } else {
58  gMessMgr->Info("##### FlowPhiWgt: FlowEvent pointer null");
59  return kStOK;
60  }
61  if (Debug()) StMaker::PrintInfo();
62 
63  return kStOK;
64 }
65 
66 //-----------------------------------------------------------------------
67 
68 Int_t StFlowPhiWgtMaker::Init() {
69 
70  // Create the files
71  TString* fileName;
72  for (int n = 1; n < nCens; n++) {
73  fileName = new TString("flowPhiWgt");
74  *fileName += n;
75  fileName->Append(".root");
76  phiWgtFile[n] = new TFile(fileName->Data(), "RECREATE");
77  delete fileName;
78  }
79 
80  // Book histograms
81  const float phiMin = 0.;
82  const float phiMax = twopi;
83  TString* histTitle;
84 
85  //ZDCSMD Phi Weight
86  mHistZDCSMDPsiWgtEast = new TH1F("Flow_ZDCSMDPsiWgtEast","Flow_ZDCSMDPsiWgtEast",
87  Flow::zdcsmd_nPsiBins,-twopi/2.,twopi/2.);
88  mHistZDCSMDPsiWgtWest = new TH1F("Flow_ZDCSMDPsiWgtWest","Flow_ZDCSMDPsiWgtWest",
89  Flow::zdcsmd_nPsiBins,-twopi/2.,twopi/2.);
90 
91  // for each centrality
92  for (int n = 1; n < nCens; n++) {
93  phiWgtFile[n]->cd(); // each file is a directory
94 
95  // for each selection
96  for (int k = 0; k < Flow::nSels; k++) {
97 
98  // for each harmonic
99  for (int j = 0; j < 2; j++) {
100 
101  // Phi lab
102  // Tpc (FarEast)
103  histTitle = new TString("Flow_Phi_FarEast_Sel");
104  *histTitle += k+1;
105  histTitle->Append("_Har");
106  *histTitle += j+1;
107  hist[k].histCen[n].histHar[j].mHistPhiFarEast =
108  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
109  hist[k].histCen[n].histHar[j].mHistPhiFarEast->SetXTitle
110  ("Azimuthal Angles (rad)");
111  hist[k].histCen[n].histHar[j].mHistPhiFarEast->SetYTitle("Counts");
112  delete histTitle;
113 
114  // Tpc (East)
115  histTitle = new TString("Flow_Phi_East_Sel");
116  *histTitle += k+1;
117  histTitle->Append("_Har");
118  *histTitle += j+1;
119  hist[k].histCen[n].histHar[j].mHistPhiEast =
120  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
121  hist[k].histCen[n].histHar[j].mHistPhiEast->SetXTitle
122  ("Azimuthal Angles (rad)");
123  hist[k].histCen[n].histHar[j].mHistPhiEast->SetYTitle("Counts");
124  delete histTitle;
125 
126  // Tpc (West)
127  histTitle = new TString("Flow_Phi_West_Sel");
128  *histTitle += k+1;
129  histTitle->Append("_Har");
130  *histTitle += j+1;
131  hist[k].histCen[n].histHar[j].mHistPhiWest =
132  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
133  hist[k].histCen[n].histHar[j].mHistPhiWest->SetXTitle
134  ("Azimuthal Angles (rad)");
135  hist[k].histCen[n].histHar[j].mHistPhiWest->SetYTitle("Counts");
136  delete histTitle;
137 
138  // Tpc (FarWest)
139  histTitle = new TString("Flow_Phi_FarWest_Sel");
140  *histTitle += k+1;
141  histTitle->Append("_Har");
142  *histTitle += j+1;
143  hist[k].histCen[n].histHar[j].mHistPhiFarWest =
144  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
145  hist[k].histCen[n].histHar[j].mHistPhiFarWest->SetXTitle
146  ("Azimuthal Angles (rad)");
147  hist[k].histCen[n].histHar[j].mHistPhiFarWest->SetYTitle("Counts");
148  delete histTitle;
149 
150  // Ftpc (FarEast)
151  histTitle = new TString("Flow_Phi_FtpcFarEast_Sel");
152  *histTitle += k+1;
153  histTitle->Append("_Har");
154  *histTitle += j+1;
155  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarEast =
156  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
157  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarEast->SetXTitle
158  ("Azimuthal Angles (rad)");
159  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarEast->SetYTitle("Counts");
160  delete histTitle;
161 
162  // Ftpc (East)
163  histTitle = new TString("Flow_Phi_FtpcEast_Sel");
164  *histTitle += k+1;
165  histTitle->Append("_Har");
166  *histTitle += j+1;
167  hist[k].histCen[n].histHar[j].mHistPhiFtpcEast =
168  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
169  hist[k].histCen[n].histHar[j].mHistPhiFtpcEast->SetXTitle
170  ("Azimuthal Angles (rad)");
171  hist[k].histCen[n].histHar[j].mHistPhiFtpcEast->SetYTitle("Counts");
172  delete histTitle;
173 
174  // Ftpc (West)
175  histTitle = new TString("Flow_Phi_FtpcWest_Sel");
176  *histTitle += k+1;
177  histTitle->Append("_Har");
178  *histTitle += j+1;
179  hist[k].histCen[n].histHar[j].mHistPhiFtpcWest =
180  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
181  hist[k].histCen[n].histHar[j].mHistPhiFtpcWest->SetXTitle
182  ("Azimuthal Angles (rad)");
183  hist[k].histCen[n].histHar[j].mHistPhiFtpcWest->SetYTitle("Counts");
184  delete histTitle;
185 
186  // Ftpc (FarWest)
187  histTitle = new TString("Flow_Phi_FtpcFarWest_Sel");
188  *histTitle += k+1;
189  histTitle->Append("_Har");
190  *histTitle += j+1;
191  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarWest =
192  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
193  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarWest->SetXTitle
194  ("Azimuthal Angles (rad)");
195  hist[k].histCen[n].histHar[j].mHistPhiFtpcFarWest->SetYTitle("Counts");
196  delete histTitle;
197 
198 
199  // PhiWgt new
200  // Tpc (FarEast)
201  histTitle = new TString("Flow_Phi_Weight_FarEast_Sel");
202  *histTitle += k+1;
203  histTitle->Append("_Har");
204  *histTitle += j+1;
205  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast =
206  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
207  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->Sumw2();
208  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->SetXTitle
209  ("Azimuthal Angles (rad)");
210  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->SetYTitle("PhiWgt");
211  delete histTitle;
212 
213  // Tpc (East)
214  histTitle = new TString("Flow_Phi_Weight_East_Sel");
215  *histTitle += k+1;
216  histTitle->Append("_Har");
217  *histTitle += j+1;
218  hist[k].histCen[n].histHar[j].mHistPhiWgtEast =
219  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
220  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->Sumw2();
221  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->SetXTitle
222  ("Azimuthal Angles (rad)");
223  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->SetYTitle("PhiWgt");
224  delete histTitle;
225 
226  // Tpc (West)
227  histTitle = new TString("Flow_Phi_Weight_West_Sel");
228  *histTitle += k+1;
229  histTitle->Append("_Har");
230  *histTitle += j+1;
231  hist[k].histCen[n].histHar[j].mHistPhiWgtWest =
232  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
233  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->Sumw2();
234  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->SetXTitle
235  ("Azimuthal Angles (rad)");
236  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->SetYTitle("PhiWgt");
237  delete histTitle;
238 
239  // Tpc (FarWest)
240  histTitle = new TString("Flow_Phi_Weight_FarWest_Sel");
241  *histTitle += k+1;
242  histTitle->Append("_Har");
243  *histTitle += j+1;
244  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest =
245  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBins, phiMin, phiMax);
246  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->Sumw2();
247  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->SetXTitle
248  ("Azimuthal Angles (rad)");
249  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->SetYTitle("PhiWgt");
250  delete histTitle;
251 
252  // Ftpc (FarEast)
253  histTitle = new TString("Flow_Phi_Weight_FtpcFarEast_Sel");
254  *histTitle += k+1;
255  histTitle->Append("_Har");
256  *histTitle += j+1;
257  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast =
258  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
259  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->Sumw2();
260  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->SetXTitle
261  ("Azimuthal Angles (rad)");
262  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->SetYTitle("PhiWgt");
263  delete histTitle;
264 
265  // Ftpc (East)
266  histTitle = new TString("Flow_Phi_Weight_FtpcEast_Sel");
267  *histTitle += k+1;
268  histTitle->Append("_Har");
269  *histTitle += j+1;
270  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast =
271  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
272  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->Sumw2();
273  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->SetXTitle
274  ("Azimuthal Angles (rad)");
275  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->SetYTitle("PhiWgt");
276  delete histTitle;
277 
278  // Ftpc (West)
279  histTitle = new TString("Flow_Phi_Weight_FtpcWest_Sel");
280  *histTitle += k+1;
281  histTitle->Append("_Har");
282  *histTitle += j+1;
283  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest =
284  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
285  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->Sumw2();
286  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->SetXTitle
287  ("Azimuthal Angles (rad)");
288  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->SetYTitle("PhiWgt");
289  delete histTitle;
290 
291  // Ftpc (FarWest)
292  histTitle = new TString("Flow_Phi_Weight_FtpcFarWest_Sel");
293  *histTitle += k+1;
294  histTitle->Append("_Har");
295  *histTitle += j+1;
296  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest =
297  new TH1D(histTitle->Data(), histTitle->Data(), Flow::nPhiBinsFtpc, phiMin, phiMax);
298  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->Sumw2();
299  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->SetXTitle
300  ("Azimuthal Angles (rad)");
301  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->SetYTitle("PhiWgt");
302  delete histTitle;
303 
304  }
305  }
306  }
307 
308  gMessMgr->SetLimit("##### FlowPhiWgt", 2);
309  gMessMgr->Info("##### FlowPhiWgt: $Id: StFlowPhiWgtMaker.cxx,v 1.7 2004/12/09 23:47:09 posk Exp $");
310 
311  return StMaker::Init();
312 }
313 
314 //-----------------------------------------------------------------------
315 
316 void StFlowPhiWgtMaker::FillParticleHistograms() {
317  // Fill histograms with event quantities
318 
319  mHistZDCSMDPsiWgtEast->Fill(pFlowEvent->ZDCSMD_PsiEst());
320  mHistZDCSMDPsiWgtWest->Fill(pFlowEvent->ZDCSMD_PsiWst());
321 
322  // Fill histograms from the particles
323 
324  // Centrality and vertex of this event
325  int iCen = pFlowEvent->Centrality();
326  if (!iCen) { return; } // skip centrality = 0
327  Float_t vertexZ = pFlowEvent->VertexPos().z();
328 
329  // Initialize Iterator
330  StFlowTrackCollection* pFlowTracks = pFlowEvent->TrackCollection();
331  StFlowTrackIterator itr;
332 
333  for (itr = pFlowTracks->begin(); itr != pFlowTracks->end(); itr++) {
334  StFlowTrack* pFlowTrack = *itr;
335 
336  float phi = pFlowTrack->Phi();
337  if (phi < 0.) phi += twopi;
338  float eta = pFlowTrack->Eta();
339  float pt = pFlowTrack->Pt();
340  float zFirstPoint = 0.;
341  float zLastPoint = 0.;
342  if (pFlowEvent->FirstLastPoints()) {
343  zFirstPoint = pFlowTrack->ZFirstPoint();
344  zLastPoint = pFlowTrack->ZLastPoint();
345  }
346  StTrackTopologyMap map = pFlowTrack->TopologyMap();
347 
348  for (int k = 0; k < Flow::nSels; k++) {
349  pFlowSelect->SetSelection(k);
350  for (int j = 0; j < 2; j++) {
351  bool oddHar = (j+1) % 2;
352  pFlowSelect->SetHarmonic(j);
353  if (pFlowSelect->Select(pFlowTrack)) {
354  // Get detID
355  StDetectorId detId;
356  Bool_t kTpcFarEast = kFALSE;
357  Bool_t kTpcEast = kFALSE;
358  Bool_t kTpcWest = kFALSE;
359  Bool_t kTpcFarWest = kFALSE;
360  Bool_t kFtpcFarEast = kFALSE;
361  Bool_t kFtpcEast = kFALSE;
362  Bool_t kFtpcWest = kFALSE;
363  Bool_t kFtpcFarWest = kFALSE;
364  if (map.hasHitInDetector(kTpcId) || (map.data(0) == 0 && map.data(1) == 0)) {
365  // Tpc track, or TopologyMap not available
366  detId = kTpcId;
367  // Set TpcEast and West
368  if (pFlowEvent->FirstLastPoints()) {
369  if (zFirstPoint > 0. && zLastPoint > 0.) {
370  kTpcFarWest = kTRUE;
371  } else if (zFirstPoint > 0. && zLastPoint < 0.) {
372  kTpcWest = kTRUE;
373  } else if (zFirstPoint < 0. && zLastPoint > 0.) {
374  kTpcEast = kTRUE;
375  } else {
376  kTpcFarEast = kTRUE;
377  }
378  } else {
379  if (eta > 0. && vertexZ > 0.) {
380  kTpcFarWest = kTRUE;
381  } else if (eta > 0. && vertexZ < 0.) {
382  kTpcWest = kTRUE;
383  } else if (eta < 0. && vertexZ > 0.) {
384  kTpcEast = kTRUE;
385  } else {
386  kTpcFarEast = kTRUE;
387  }
388  }
389  } else if (map.trackFtpcEast()) {
390  detId = kFtpcEastId; // eta < 0.
391  if (vertexZ > 0.) {
392  kFtpcEast = kTRUE;
393  } else { // vertexZ < 0.
394  kFtpcFarEast = kTRUE;
395  }
396  } else if (map.trackFtpcWest()) {
397  detId = kFtpcWestId; // eta > 0.
398  if (vertexZ > 0.) {
399  kFtpcFarWest = kTRUE;
400  } else { // vertexZ < 0.
401  kFtpcWest = kTRUE;
402  }
403  } else {
404  detId = kUnknownId;
405  }
406 
407  // Calculate weights for filling histograms
408  float wt = 1.;
409  if (pFlowEvent->PtWgt()) {
410  wt *= (pt < pFlowEvent->PtWgtSaturation()) ? pt : pFlowEvent->PtWgtSaturation(); // pt weighting going constant
411  }
412  float etaAbs = fabs(eta);
413  if (pFlowEvent->EtaWgt() && oddHar && etaAbs > 1.) { wt *= etaAbs; }
414 
415  // Fill histograms
416  if (kFtpcFarEast) {
417  hist[k].histCen[iCen].histHar[j].mHistPhiFtpcFarEast->Fill(phi,wt);
418  } else if (kFtpcEast) {
419  hist[k].histCen[iCen].histHar[j].mHistPhiFtpcEast->Fill(phi,wt);
420  } else if (kFtpcWest) {
421  hist[k].histCen[iCen].histHar[j].mHistPhiFtpcWest->Fill(phi,wt);
422  } else if (kFtpcFarWest) {
423  hist[k].histCen[iCen].histHar[j].mHistPhiFtpcFarWest->Fill(phi,wt);
424  } else if (kTpcFarEast){
425  hist[k].histCen[iCen].histHar[j].mHistPhiFarEast->Fill(phi,wt);
426  } else if (kTpcEast){
427  hist[k].histCen[iCen].histHar[j].mHistPhiEast->Fill(phi,wt);
428  } else if (kTpcWest){
429  hist[k].histCen[iCen].histHar[j].mHistPhiWest->Fill(phi,wt);
430  } else if (kTpcFarWest){
431  hist[k].histCen[iCen].histHar[j].mHistPhiFarWest->Fill(phi,wt);
432  }
433 
434  }
435  }
436  }
437  }
438 
439 }
440 
441 //-----------------------------------------------------------------------
442 
444  // Outputs phiWgt values
445 
446  // PhiWgt histogram collections
447  TOrdCollection* phiWgtHistNames[nCens];
448  for (int n = 1; n < nCens; n++) {
449  phiWgtHistNames[n] = new TOrdCollection(Flow::nSels*2 + 2);
450  }
451 
452  cout << endl << "##### PhiWgt Maker:" << endl;
453 
454  for (int n = 1; n < nCens; n++) {
455  for (int k = 0; k < Flow::nSels; k++) {
456  for (int j = 0; j < 2; j++) {
457  // Calculate PhiWgt
458  double meanFarEast = hist[k].histCen[n].histHar[j].mHistPhiFarEast->
459  Integral() / (double)Flow::nPhiBins;
460  double meanEast = hist[k].histCen[n].histHar[j].mHistPhiEast->
461  Integral() / (double)Flow::nPhiBins;
462  double meanWest = hist[k].histCen[n].histHar[j].mHistPhiWest->
463  Integral() / (double)Flow::nPhiBins;
464  double meanFarWest = hist[k].histCen[n].histHar[j].mHistPhiFarWest->
465  Integral() / (double)Flow::nPhiBins;
466  double meanFtpcFarEast = hist[k].histCen[n].histHar[j].mHistPhiFtpcFarEast->
467  Integral() / (double)Flow::nPhiBinsFtpc;
468  double meanFtpcEast = hist[k].histCen[n].histHar[j].mHistPhiFtpcEast->
469  Integral() / (double)Flow::nPhiBinsFtpc;
470  double meanFtpcWest = hist[k].histCen[n].histHar[j].mHistPhiFtpcWest->
471  Integral() / (double)Flow::nPhiBinsFtpc;
472  double meanFtpcFarWest = hist[k].histCen[n].histHar[j].mHistPhiFtpcFarWest->
473  Integral() / (double)Flow::nPhiBinsFtpc;
474 
475  // Tpc
476  for (int i = 0; i < Flow::nPhiBins; i++) {
477  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->
478  SetBinContent(i+1,meanFarEast);
479  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->
480  SetBinError(i+1, 0.);
481  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->
482  SetBinContent(i+1, meanEast);
483  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->
484  SetBinError(i+1, 0.);
485  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->
486  SetBinContent(i+1, meanWest);
487  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->
488  SetBinError(i+1, 0.);
489  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->
490  SetBinContent(i+1,meanFarWest);
491  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->
492  SetBinError(i+1, 0.);
493  }
494 
495  // Ftpc
496  for (int i = 0; i < Flow::nPhiBinsFtpc; i++) {
497  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->
498  SetBinContent(i+1, meanFtpcFarEast);
499  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->
500  SetBinError(i+1, 0.);
501  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->
502  SetBinContent(i+1, meanFtpcEast);
503  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->
504  SetBinError(i+1, 0.);
505  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->
506  SetBinContent(i+1, meanFtpcWest);
507  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->
508  SetBinError(i+1, 0.);
509  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->
510  SetBinContent(i+1, meanFtpcFarWest);
511  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->
512  SetBinError(i+1, 0.);
513  }
514 
515  // Tpc
516  hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast->
517  Divide(hist[k].histCen[n].histHar[j].mHistPhiFarEast);
518  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFarEast);
519  hist[k].histCen[n].histHar[j].mHistPhiWgtEast->
520  Divide(hist[k].histCen[n].histHar[j].mHistPhiEast);
521  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtEast);
522  hist[k].histCen[n].histHar[j].mHistPhiWgtWest->
523  Divide(hist[k].histCen[n].histHar[j].mHistPhiWest);
524  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtWest);
525  hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest->
526  Divide(hist[k].histCen[n].histHar[j].mHistPhiFarWest);
527  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFarWest);
528 
529  // Ftpc
530  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast->
531  Divide(hist[k].histCen[n].histHar[j].mHistPhiFtpcFarEast);
532  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarEast);
533  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast->
534  Divide(hist[k].histCen[n].histHar[j].mHistPhiFtpcEast);
535  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcEast);
536  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest->
537  Divide(hist[k].histCen[n].histHar[j].mHistPhiFtpcWest);
538  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcWest);
539  hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest->
540  Divide(hist[k].histCen[n].histHar[j].mHistPhiFtpcFarWest);
541  phiWgtHistNames[n]->AddLast(hist[k].histCen[n].histHar[j].mHistPhiWgtFtpcFarWest);
542  }
543  }
544  phiWgtHistNames[n]->AddLast(mHistZDCSMDPsiWgtEast);
545  phiWgtHistNames[n]->AddLast(mHistZDCSMDPsiWgtWest);
546  }
547  //GetHistList()->ls();
548 
549  // Make text object
550  TText* textInfo = 0;
551  if (pFlowEvent->FirstLastPoints()) {
552  char chInfo[400];
553  sprintf(chInfo, "%s%d%s%d%s", " pt weight= ", pFlowEvent->PtWgt(),
554  ", eta weight= ", pFlowEvent->EtaWgt(), "\n");
555  textInfo = new TText(0,0,chInfo);
556  }
557 
558  // Write PhiWgt histograms
559  for (int n = 1; n < nCens; n++) {
560  phiWgtFile[n]->cd();
561  if (pFlowEvent->FirstLastPoints()) { textInfo->Write("info"); }
562  phiWgtHistNames[n]->Write();
563  phiWgtFile[n]->Close();
564  delete phiWgtHistNames[n];
565  }
566 
567  if (pFlowEvent->FirstLastPoints()) delete textInfo;
568 
569  delete pFlowSelect;
570 
571  return StMaker::Finish();
572 }
573 
574 //-----------------------------------------------------------------------
575 
577 //
578 // $Log: StFlowPhiWgtMaker.cxx,v $
579 // Revision 1.7 2004/12/09 23:47:09 posk
580 // Minor changes in code formatting.
581 // Added hist for TPC primary dca to AnalysisMaker.
582 //
583 // Revision 1.6 2004/12/07 23:10:22 posk
584 // Only odd and even phiWgt hists. If the old phiWgt file contains more than
585 // two harmonics, only the first two are read. Now writes only the first two.
586 //
587 // Revision 1.5 2004/08/24 20:22:40 oldi
588 // Minor modifications to avoid compiler warnings.
589 //
590 // Revision 1.4 2004/05/31 20:09:25 oldi
591 // PicoDst format changed (Version 7) to hold ZDC SMD information.
592 // Trigger cut modified to comply with TriggerCollections.
593 // Centrality definition for 62 GeV data introduced.
594 // Minor bug fixes.
595 //
596 // Revision 1.3 2003/09/02 17:58:11 perev
597 // gcc 3.2 updates + WarnOff
598 //
599 // Revision 1.2 2003/07/30 22:01:55 oldi
600 // PtWgtSaturation parameter introduced.
601 //
602 // Revision 1.1 2003/05/16 20:44:49 posk
603 // First commit of StFlowPhiWgtMaker
604 //
605 //
Definition: Stypes.h:40
virtual Int_t Finish()
Definition: StMaker.cxx:776