HYDRA_development_version
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
hparticleevtchara.cc
Go to the documentation of this file.
1 //_HADES_CLASS_DESCRIPTION ////////////////////////////////////////////////////////////////////////////
2 //*-- AUTHOR : B. Kardan 28.08.2019
3 //*-- VERSION : 0.73
4 //
5 // HParticleEvtChara
6 //
7 // Purpose: EventCharacterization
8 // - Centrality from Hit(TOF,RPC,FW) and Track Estimators in Data/Sim
9 // - QVector and Phi (Reaction Plane Estimate) from FW
10 // - Event-weight for downscaled Events PT2
11 //
12 //--------------------------------------------------------------------------
13 // Usage:
14 //
15 // - input files can be found at : example /cvmfs/hades.gsi.de/param/eventchara/
16 //
17 // - to define the ParameterFile where Classes and Estimators are stored use
18 // setParameterFile("/cvmfs/hades.gsi.de/param/eventchara/ParameterFile.root")
19 //
20 // - to print the definition of estimator & class use
21 // printCentralityClass(HParticleEvtChara::kTOFRPC, HParticleEvtChara::k10);
22 //
23 // - to get the CentralityClass of an event (with estimator and class definition) use
24 // getCentralityClass(HParticleEvtChara::kTOFRPC, HParticleEvtChara::k10);
25 //
26 // - to get the CentralityPercentile of an event (only estimator is needed) use
27 // getCentralityPercentile(HParticleEvtChara::kTOFRPC);
28 //
29 // - to get the EventWeight to re-weight downscaled events use
30 // getEventWeight();
31 //
32 // - to get the EventPlane with ReCentering use
33 // getEventPlane(HParticleEvtChara::kDefault);
34 //
35 //--------------------------------------------------------------------------
36 // Estimators:
37 // TOFRPC - (default) TOF and RPC hit multiplicity in time-window
38 // TOF - TOF hit multiplicity in time-window
39 // RPC - RPC hit multiplicity in time-window
40 // TOFRPCtot - total TOF and RPC hit multiplicity in event-window
41 // SelectedParticleCand - selected Particle multiplicity
42 // PrimaryParticleCand - primary Particle multiplicity
43 // SelectedParticleCandCorr - selected Particle multiplicity corrected by the
44 // running mean and scaled with the reference mean
45 // (selTrack * referenceMean/<selTrack>)
46 // SelectedParticleCandSecCorr
47 // - selected Particle multiplicity corrected by the
48 // running mean and scaled with the reference mean
49 // (selTrack * referenceMean/<selTrack>)
50 // SelectedParticleCandNorm - selected Particle multiplicity normalized by
51 // the running mean (selTrack/<selTrack>)
52 // FWSumChargeSpec - sum of charge (dE/dx in a.u.) of Spectator in
53 // FW acceptance with beta~0.9
54 // FWSumChargeZ - sum of charge (int Z up till charge-state 14
55 // with individuel fixed cuts in dE/dx each FW-cell)
56 // of Spectator in FW acceptance with beta~0.9
57 //
58 // Classes:
59 // 2 - 2% classes
60 // 5 - 5% classes
61 // 10 - (default) 10% classes
62 // 13 - 13% classes
63 // 15 - 15% classes
64 // 20 - 20% classes
65 // FOPI - FOPI centrality classes
66 //
67 // EventPlaneCorrection:
68 // kNoCorrection - EP only selection of spectator candidates in FW,
69 // no Correction
70 // kSmear - smearing of x and y of each FW hit inside cell size
71 // kShiftFW - global shift x and y with the centre of gravity and
72 // scale x and y with the sigma of the distribution of
73 // all FW hits in all events (per class/day)
74 // kWeightCharge - weigthing with charge-state up to 14 with individuel
75 // fixed timing- and dE/dx-cuts each FW-cell
76 // kReCentering - re-centering of QVector with <Qx> and <Qy> (calc. evt-by-evt)
77 // (only first harmonic correction)
78 // kScaling - scaling of QVector with the sigma of Qx and Qy (calc. evt-by-evt)
79 // kRotation - rotation of EP via residual Fourier harmonics up to 8 cos and sin terms
80 // after FWshift and FWscaling
81 //
82 // kDefault - recommended option (kShiftFW|kWeightCharge|kRotation)
83 //--------------------------------------------------------------------------
84 // quick how-To:
85 // - for an example see demo-macro:
86 // /scripts/batch/GE/loopDST/loopDST.C
87 //
88 // 1. include header-file
89 // #include "HParticleEvtChara.h"
90 //
91 // 2. if you use Hloop check that following catagories are loaded:
92 //
93 // HParticleEvtInfo
94 // HParticleCand
95 // HWallHit
96 //
97 // if(!loop.setInput("-*,+HParticleEvtInfo, +HParticleCand, +HWallHit")){
98 // cout<<"READBACK: ERROR : cannot read input !"<<endl;
99 // }
100 //
101 // 3. before eventLoop:
102 //
103 // HParticleEvtChara evtChara;
104 // TString ParameterfileCVMFS = "/lustre/nyx/hades/user/bkardan/param/centrality_epcorr_apr12_gen8_2019_02_pass30.root";
105 //
106 // //due to the overlap of day126 there is dedidacate param-file for the reverse-field runs
107 // if(isRevFieldDATA) ParameterfileCVMFS = "/lustre/nyx/hades/user/bkardan/param/centrality_epcorr_apr12_gen8_revfield_2019_02_pass29.root";
108 //
109 // //if you need centrality-values for sim UrQMD gen8a you find them here
110 // //due missing primary light nuclei in UrQMD... there is no correction on the FW-EventPlane in this version!
111 // if(isSimulation) ParameterfileCVMFS = "/lustre/nyx/hades/user/bkardan/param/centrality_epcorr_sim_au1230au_gen8a_UrQMD_minbias_2019_04_pass0.root";
112 // //if(isSimulation) ParameterfileCVMFS = "/lustre/nyx/hades/user/bkardan/param/centrality_epcorr_sim_au1230au_gen9vertex_UrQMD_minbias_2019_04_pass0.root";
113 //
114 // if(!evtChara.setParameterFile(ParameterfileCVMFS)){
115 // cout << "Parameterfile not found !!! " << endl;return kFALSE;
116 // }
117 //
118 // if(!evtChara.init()) {
119 // cout << "HParticleEvtChara not init!!! " << endl;return kFALSE;
120 // }
121 //
122 // Int_t eCentEst = HParticleEvtChara::kTOFRPC;
123 // Int_t eCentClass = HParticleEvtChara::k10;
124 // Int_t eEPcorr = HParticleEvtChara::kDefault;
125 // cout << "\t selected EPcorrection method is: " << evtChara.getStringEventPlaneCorrection(eEPcorr) << endl;
126 //
127 //
128 // evtChara.printCentralityClass(eCentEst, eCentClass);
129 //
130 // 4. inside event-loop:
131 //
132 // Float_t event_weight = evtChara.getEventWeight(); // event_weight dependent if PT2(down-scaled) or PT3
133 //
134 // Int_t CentralityClassTOFRPC = evtChara.getCentralityClass(eCentEst, eCentClass);
135 // // 10% Centrality-Classes: 1(0-10%) - 5(40-50%) ... 0:Overflow max:Underflow
136 //
137 // Float_t CentralityTOFRPC = evtChara.getCentralityPercentile(eCentClass);
138 // // CentralityPercentile: 0 - 100% percentile of the total cross section
139 // // 101% Over-,Underflow or Outlier Events
140 //
141 // Float_t EventPlane = evtChara.getEventPlane(eEPcorr);
142 // //EventPlane: 0 - 2pi (in rad) re-centered & scaled EP
143 // // -1 no EP could be reconstructed
144 //
145 // //for EP-resolution use the sub-event:
146 // Float_t EventPlaneA = evtChara.getEventPlane(eEPcorr,1);
147 // Float_t EventPlaneB = evtChara.getEventPlane(eEPcorr,2);
148 //
149 // //check if the EP(and subEvent EP) could be reconstructed, most of the case less than 4 Hits in FW
150 // if(EventPlane == -1) continue;
151 // if(EventPlaneA == -1 || EventPlaneB == -1 ) continue;
152 //
153 // check if you have any further selections on FW-hit multiplicites in you own code,
154 // this will bias the result!
155 //
156 // 5. to calculate the EP-resolution one option is the following:
157 // (Ollitrault method. - approximation valid for high chi-values)
158 // inside event-loop fill the difference of the two sub-event EP
159 // into a histogramm for each centrality class:
160 //
161 // Float_t deltaPhiAB = TMath::Abs(TVector2::Phi_mpi_pi(EventPlaneA-EventPlaneB));
162 // h->Fill(deltaPhiAB);
163 //
164 // after your event-analysis you can finalize your results by calculating the EP-resolution for each centrality class
165 // by the ratio of event N(pi/2 < deltaEP_AB < pi)/N(0 < deltaEP_AB < pi):
166 //
167 // Int_t i90 = h->FindBin(0.5*TMath::Pi());
168 // Int_t i180 = h->FindBin(TMath::Pi());
169 // Double_t ratio = h->Integral(i90,i180) /h->Integral(1,i180);
170 // Double_t dChi = sqrt(-2.*TMath::Log(2.*ratio));
171 // Double_t dEPres1 = HParticleEvtChara::ComputeResolution( dChi ); // for first harmonic
172 // Double_t dEPres2 = HParticleEvtChara::ComputeResolution( dChi , 2); // for second ...
173 //
174 // 6. the second option will be described here soon here ;)
175 // (Voloshin method - precise iterative resolution search)
176 // inside event-loop fill the Cosine of the difference of the two sub-event EP
177 // into a histogramm for each centrality class:
178 //
179 // Float_t CosDeltaPhiAB = TMath::Cos(EventPlaneA-EventPlaneB);
180 // h2->Fill(CosDeltaPhiAB);
181 //
182 // after your event-analysis you can finalize your results by calculating the EP-resolution for each centrality class
183 // by the mean-value:
184 //
185 // Double_t MeanCosDeltaPhiAB = h2->GetMean(); //< cos n delta_phi>;
186 // Double_t dV = TMath::Sqrt(MeanCosDeltaPhiAB);
187 // Double_t dChi = FindXi(dV,1e-6, 1); // for first harmonic
188 // Double_t dEPres1 = HParticleEvtChara::ComputeResolution( dChi ); // for first harmonic
189 // Double_t dEPres2 = HParticleEvtChara::ComputeResolution( dChi , 2); // for second ...
190 // dV = ComputeResolution( TMath::Sqrt2()*dChi , 1);
191 // printf("An estimate of the event plane resolution first order is: %f\n", dV );
192 //
193 //
194 // 7. to correct your flow-values use 1./dEPres1
195 //
196 // 8. test your macro with data and if you real have the right parameterfile
197 // with right values for centrality and EPcorections, check the output of
198 //
199 // evtChara.printCentralityClass(eCentEst, eCentClass);
200 //
201 // 9. and now happy plotting ;)
202 //
203 //--------------------------------------------------------------------------
204 // Change History:
205 //
206 // 19.02.2019 release of 0.7
207 // 24.04.2019 implementation of useFWCut() and check if FWcut hist are loaded
208 // 04.06.2019 fix of embeded-sim into data and removed of high-momentum condition
209 ////////////////////////////////////////////////////////////////////////////
210 
211 #include "hparticleevtchara.h"
213 
214 
215 HParticleEvtChara::HParticleEvtChara(const Text_t* name,const Text_t* title)
216  : HReconstructor(name,title) {
217  fParticleEvtInfoCat = 0;
218  fParameterFile = "centrality.root";
219  fPostFix = "";
220  isSimulation = kFALSE;
221  fReferenceMeanSelTrack= 36.44;
222  fDayOfYear = 0;
223  fEventPlaneCorrectionFlag = kNoCorrection;
224  fQVectorCalcDone = kFALSE;
225  //maxFWCells = 302;
226 
227  useFWCut.resize(kNumFWCutValues);
228  useFWCut[0]=kFALSE; //kBetaCuts
229  useFWCut[1]=kTRUE; //kTimeCuts
230  useFWCut[2]=kTRUE; //kChargeCuts
231  excludeNoisyFWcells = kTRUE;
232 
233  fRandom = new TRandom2();
234 
235  vQPhi.resize(3);
236 
237  fCentralityPercentileHist.resize(kNumCentralityEstimator);
238  fCentralityHist.resize(kNumCentralityEstimator);
239  fEstimatorHist.resize(kNumCentralityEstimator);
240  for (int centEst = 0; centEst < (int)kNumCentralityEstimator; ++centEst) fCentralityHist[centEst].resize(kNumCentralityClass);
241  fEventPlaneCorrectionHist.resize(kNumEventPlaneParameter);
242  fFWCutValuesHist.resize(kNumFWCutValues);
243  for (int cutValue = 0; cutValue < (int)kNumFWCutValues; ++cutValue) fFWCutValuesHist[cutValue].resize(MAXFWCELLS);
244 
245 
246  fFWminBeta.resize(3);
247  fFWmaxBeta.resize(3);
248  fFWminCharge.resize(3);
249  fFWminBeta[0]=0.84; fFWmaxBeta[0]=1.2; fFWminCharge[0]=80; // small cells
250  fFWminBeta[1]=0.85; fFWmaxBeta[1]=1.2; fFWminCharge[1]=84; // medium cells
251  fFWminBeta[2]=0.8; fFWmaxBeta[2]=1.2; fFWminCharge[2]=88; // large cells
252 
253  fFWChargeCuts.resize(6);
254  fFWChargeCuts[5]=386;
255  fFWChargeCuts[4]=339;
256  fFWChargeCuts[3]=298;
257  fFWChargeCuts[2]=241;
258  fFWChargeCuts[1]=175;
259  fFWChargeCuts[0]=84 ;
260 
261  }
262 
264 }
265 
267 {
268  if(gHades){
270  if(evt){
271 
272  HCategory* catKin=HCategoryManager::getCategory(catGeantKine,kFALSE,"catGeantKine");
273  if(catKin) {isSimulation=kTRUE; Info("init()","GeantKine found - is Simulation");}
274  else {isSimulation=kFALSE;Info("init()","GeantKine not found - is not Simulation");}
275 
276 
277  fParticleEvtInfoCat = HCategoryManager::getCategory(catParticleEvtInfo,kTRUE,"catParticleEvtInfo, from HParticleEvtChara::init()");
278  if(!fParticleEvtInfoCat) { Info("init()","No catParticleEvtInfo in input!"); return kFALSE;}
279 
280  fParticleCandCat = HCategoryManager::getCategory(catParticleCand,kTRUE,"catParticleCand, from HParticleEvtChara::init()");
281  if(!fParticleCandCat) { Info("init()","No catParticleCand in input!");}
282 
283  fCatWallHit = HCategoryManager::getCategory(catWallHit,kTRUE,"catWallHit, from HParticleEvtChara::init()");
284  if(!fCatWallHit) { Info("init()","No catWallHit in input!");}
285 
286  } else {
287  Error("init()","NO EventStructure found!");
288  return kFALSE;
289  }
290 
291  } else {
292  Error("init()","NO HADES found!");
293  return kFALSE;
294  }
295  // read parameter file
296  if(fParameterFile){
297  if(!loadParameterFile()) return kFALSE;
298  } else {
299  Error("init()","NO Parameterfile found!");
300  return kFALSE;
301  }
302  loadDayOfYear();
303 
304  return kTRUE;
305 }
306 
308 {
309  loadDayOfYear();
310  return kTRUE;
311 }
312 
314 {
315  TString tempFileName = "";
316  if(gLoop){
317  gLoop->isNewFile(tempFileName);
318  }
319  else if(gHades){
320  if(gHades->getRuntimeDb()){
321  tempFileName = gHades->getRuntimeDb()->getCurrentFileName();
322  }
323  else if(gHades->getDataSource()){
324  tempFileName = gHades->getDataSource()->getCurrentFileName();
325  }
326  }
327  if(tempFileName==""){
328  //Error("loadDayOfYear()","File/Day not found! using default day108");
329  fDayOfYear = 108; // default FIXME in apr12
330  if(isSimulation) fDayOfYear = 0;
331  return fDayOfYear;
332  }
333  fDayOfYear = HTime::getDayFileName(HTime::stripFileName(tempFileName,kTRUE,kTRUE));
334  //Info("loadDayOfYear()",">>> Day of Year: %d",fDayOfYear);
335 
336  if( fDayOfYear>365 && isSimulation){
337  fDayOfYear = 0;
338  }
339 
340  return fDayOfYear;
341 }
342 
343 Bool_t HParticleEvtChara::setParameterFile(TString ParameterFile)
344 {
345  fParameterFile = ParameterFile;
346  TString path = gSystem->ExpandPathName(fParameterFile.Data());
347  if (gSystem->AccessPathName(path)) {
348  Error("init()","File %s does not exist!",path.Data());
349  return kFALSE;
350  }
351  return kTRUE;
352 }
353 
355 {
356  // read parameter file
357  TString path = gSystem->ExpandPathName(fParameterFile.Data());
358  if (gSystem->AccessPathName(path)) {
359  Error("loadParameterFile()","File %s does not exist!",path.Data());
360  return kFALSE;
361  } else {
362  fFile = new TFile(path,"OPEN");
363  TObject* ParameterFileVersion = (TObject*) fFile->Get("HParticleEvtCharaVersion");
364  if(ParameterFileVersion){
365  TString Version(ParameterFileVersion->GetTitle());
366  Float_t fVersion = Version.Atof();
367  if(fVersion < getVersion()){
368  Error("loadParameterFile()","File %s is out-dated! Needed Version: %02.1f -- in file: %02.1f",path.Data(), getVersion(), fVersion);
369  return kFALSE;
370  }
371  cout<<"\n--------------------------------------------------------------------------------------" << endl;
372  Info("loadParameterFile()",">>> Parameter input file (ver. %02.1f) : %s",fVersion,path.Data());
373 
374  cout<<"\n--------------------------------------------------------------------------------------" << endl;
375  }
376  else{
377  Error("loadParameterFile()","In File %s no Version information found!",path.Data());
378  return kFALSE;
379  }
380  }
384  return kTRUE;
385 }
386 
388 {
389  // read parameter file
390  TString path = gSystem->ExpandPathName(fParameterFile.Data());
391  //if (gSystem->AccessPathName(path)) {
392  // Error("loadParameterFile()","File %s does not exist!",path.Data());
393  // return kFALSE;
394  //} else {
395  fFile = new TFile(path,"CREATE");
396  if(!fFile) return kFALSE;
397  cout << "Version of HParticleEvtChara: " << getVersion() << endl;
398  TString version=Form("%02.1f",getVersion());
399  TNamed OutputVersion("HParticleEvtCharaVersion",version.Data());
400  OutputVersion.Write();
404  return kTRUE;
405 }
406 
408 {
409 
410  cout<<"\n--------------------------------------------------------------------------------------" << endl;
411  Info("loadCentralityEstimatorHist()","Calibration for Centrality Estimator and Classes loading:");
412  Int_t n=0,m=0;
413  for (int centEst = 0; centEst < (int) kNumCentralityEstimator; ++centEst){
414  printf("\n>>>> %31s %s: ",getStringCentralityEstimator(centEst).Data() ,fPostFix.Data());
415  TString temp;
416  if(fPostFix.CompareTo("")==0) temp = Form("%s_percentile",getStringCentralityEstimator(centEst).Data() );
417  else temp = Form("%s_%s_percentile",getStringCentralityEstimator(centEst).Data(), fPostFix.Data() );
418  fCentralityPercentileHist[centEst] = (TH1F*) fFile->FindObjectAny( temp.Data() );
419  if(fCentralityPercentileHist[centEst]){
420  printf("percentile");
421  m++;
422  }
423 
424  for (int centC = 0; centC < (int) kNumCentralityClass; ++centC){
425  TString temp2;
426  if(fPostFix.CompareTo("")==0) temp2 = Form("%s_%s_fixedCuts",getStringCentralityEstimator(centEst).Data(), getStringCentralityClass(centC).Data() );
427  else temp2 = Form("%s_%s_%s_fixedCuts",getStringCentralityEstimator(centEst).Data(), fPostFix.Data(), getStringCentralityClass(centC).Data() );
428  fCentralityHist[centEst][centC] = (TH1F*) fFile->FindObjectAny( temp2.Data() );
429  if(fCentralityHist[centEst][centC]){
430  printf(" %5s",getStringCentralityClass(centC).Data() );
431  n++;
432  }
433  }
434  }
435  cout<< "\n\nCentrality Percentile Curves #"<< m << endl;
436  cout<< "Centrality Estimator and Classes #"<< n << endl;
437  cout<< "--------------------------------------------------------------------------------------" << endl;
438 
439  return 0;
440 }
441 
443 {
444  if(!fFile) return kFALSE;
445  fFile->mkdir("EstimatorHist");
446  fFile->cd("/EstimatorHist/");
447 
448  cout<<"\n--------------------------------------------------------------------------------------" << endl;
449  Info("saveCentralityEstimatorHist()","Estimator Hist saving:");
450  Int_t m=0;
451  for (int centEst = 0; centEst < (int) kNumCentralityEstimator; ++centEst){
452  if(fEstimatorHist[centEst]){
453  fEstimatorHist[centEst]->Write();
454  m++;
455  }
456  }
457  cout<< "\n\nEstimator Hist #"<< m << "saved" << endl;
458  cout<< "--------------------------------------------------------------------------------------" << endl;
459 
460  fFile->mkdir("Centrality");
461  fFile->cd("/Centrality/");
462 
463  cout<<"\n--------------------------------------------------------------------------------------" << endl;
464  Info("saveCentralityEstimatorHist()","Calibration for Centrality Estimator and Classes saving:");
465  Int_t n=0;
466  m=0;
467  for (int centEst = 0; centEst < (int) kNumCentralityEstimator; ++centEst){
468  printf("\n>>>> %31s %s: ",getStringCentralityEstimator(centEst).Data() ,fPostFix.Data());
469  if(fCentralityPercentileHist[centEst]){
470  fCentralityPercentileHist[centEst]->Write();
471  printf("percentile");
472  m++;
473  }
474 
475  for (int centC = 0; centC < (int) kNumCentralityClass; ++centC){
476  if(fCentralityHist[centEst][centC]){
477  fCentralityHist[centEst][centC]->Write();
478  printf(" %5s",getStringCentralityClass(centC).Data() );
479  n++;
480  }
481  }
482  }
483  cout<< "\n\nCentrality Percentile Curves #"<< m << "saved" << endl;
484  cout<< "Centrality Estimator and Classes #"<< n << "saved" << endl;
485  cout<< "--------------------------------------------------------------------------------------" << endl;
486  return 0;
487 }
488 
489 
491 {
492 
493  cout<<"\n--------------------------------------------------------------------------------------" << endl;
494  Info("loadEventPlaneCorrectionHist()","Calibration for EventPlane Correction Histogramms loading:");
495  Int_t n=0;
496  for (int epParam = 0; epParam < (int) kNumEventPlaneParameter; ++epParam){
497  printf("\n>>>> %27s : ",getStringEventPlaneParameter(epParam).Data() );
498  TString temp;
499  temp = Form("EPcorr_%s_Day_Centrality",getStringEventPlaneParameter(epParam).Data() );
500  fEventPlaneCorrectionHist[epParam] = (TProfile2D*) fFile->FindObjectAny( temp.Data() ); //FIXME
501  if(fEventPlaneCorrectionHist[epParam]){
502  //for (int dim = 0; dim < fEventPlaneCorrectionHist[qAxis]->GetDimension(); ++dim){
503  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetXaxis()->GetTitle(),
504  fEventPlaneCorrectionHist[epParam]->GetNbinsX());
505  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetYaxis()->GetTitle(),
506  fEventPlaneCorrectionHist[epParam]->GetNbinsY());
507  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetZaxis()->GetTitle(),
508  fEventPlaneCorrectionHist[epParam]->GetNbinsZ());
509  n++;
510  }
511  }
512 
513  cout<< "\n\nEventPlane Correction Histogramms #"<< n << endl;
514  cout<< "--------------------------------------------------------------------------------------" << endl;
515 
516  return 0;
517 }
518 
519 Bool_t HParticleEvtChara::addEventPlaneCorrectionHist(TProfile2D *hist, UInt_t epParam)
520 {
521  if(!hist)return kFALSE;
522  if(epParam>=kNumEventPlaneParameter) return kFALSE;
523  if(epParam==kChi) return kTRUE;
524 
525  TString temp;
526  temp = Form("EPcorr_%s_Day_Centrality",getStringEventPlaneParameter(epParam).Data() );
527  hist->SetNameTitle(temp.Data(),temp.Data());
528 
529  if(epParam==kResolution){
531 
532  temp = Form("EPcorr_%s_Day_Centrality",getStringEventPlaneParameter(kChi).Data() );
533  hist->SetNameTitle(temp.Data(),temp.Data());
534  fEventPlaneCorrectionHist[kChi] = (TH2D*) makeEPresolution(hist, kTRUE);
535  }
536  else{
537  fEventPlaneCorrectionHist[epParam] = (TProfile2D*) hist;
538  }
539  return kTRUE;
540 }
541 
543 {
544  if(!fFile) return kFALSE;
545  fFile->mkdir("EPcorr");
546  fFile->cd("/EPcorr/");
547 
548  cout<<"\n--------------------------------------------------------------------------------------" << endl;
549  Info("saveEventPlaneCorrectionHist()","Calibration for EventPlane Correction Histogramms saving:");
550  Int_t n=0;
551  for (int epParam = 0; epParam < (int) kNumEventPlaneParameter; ++epParam){
552  printf("\n>>>> %27s : ",getStringEventPlaneParameter(epParam).Data() );
553  if(fEventPlaneCorrectionHist[epParam]){
554  fEventPlaneCorrectionHist[epParam]->Write();
555  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetXaxis()->GetTitle(),
556  fEventPlaneCorrectionHist[epParam]->GetNbinsX());
557  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetYaxis()->GetTitle(),
558  fEventPlaneCorrectionHist[epParam]->GetNbinsY());
559  printf("%4s(%2d bins) ",fEventPlaneCorrectionHist[epParam]->GetZaxis()->GetTitle(),
560  fEventPlaneCorrectionHist[epParam]->GetNbinsZ());
561  n++;
562  }
563  }
564 
565  cout<< "\n\nEventPlane Correction Histogramms #"<< n << "saved" << endl;
566  cout<< "--------------------------------------------------------------------------------------" << endl;
567 
568  return 0;
569 }
570 
572 {
573  cout<<"\n--------------------------------------------------------------------------------------" << endl;
574  Info("loadFWCutValuesHist()","Histogramms with CutValues for FW Hits loading:");
575  Int_t n=0;
576  for (int ep = 0; ep < (int) kNumFWCutValues; ++ep){
577  Int_t m = 0;
578  printf("\n>>>> %27s : ",getStringFWCutValues(ep).Data() );
579  for (int iCell = 0; iCell < MAXFWCELLS ; ++iCell){
580  TString temp = Form("FWCuts_%s_cell%d",getStringFWCutValues(ep).Data(), iCell );
581  fFWCutValuesHist[ep][iCell] = (TH1F*) fFile->FindObjectAny( temp.Data() ); //FIXME
582  if(fFWCutValuesHist[ep][iCell]){
583  if(m==0){
584  printf("%4s(%2d bins) ",fFWCutValuesHist[ep][iCell]->GetXaxis()->GetTitle(),
585  fFWCutValuesHist[ep][iCell]->GetNbinsX());
586  printf("%4s(%2d bins) ",fFWCutValuesHist[ep][iCell]->GetYaxis()->GetTitle(),
587  fFWCutValuesHist[ep][iCell]->GetNbinsY());
588  }
589  n++; m++;
590  }
591  }
592  printf("active cells:%2d ", m);
593  if(m==0){
594  useFWCut[ep] = kFALSE;
595  printf(" not used for FWCut!");
596  }
597  }
598 
599  cout<< "\n\nHistogramms with CutValues for FW Hits #"<< n << endl;
600  cout<< "--------------------------------------------------------------------------------------" << endl;
601 
602  return 0;
603 }
604 
605 Bool_t HParticleEvtChara::addFWCutValuesHist(TH1 *hist, Int_t cell, UInt_t eFWCut)
606 {
607  if(!hist)return kFALSE;
608  if(eFWCut>=kNumFWCutValues) return kFALSE;
609  if(cell<0 || cell >= MAXFWCELLS) return kFALSE;
610 
611  TString temp = Form("FWCuts_%s_cell%d",getStringFWCutValues(eFWCut).Data(), cell );
612  hist->SetNameTitle(temp.Data(),temp.Data());
613  fFWCutValuesHist[eFWCut][cell] = (TH1F*) hist;
614  return kTRUE;
615 }
616 
618 {
619  if(!fFile) return kFALSE;
620  fFile->mkdir("FWCutValue");
621  fFile->cd("/FWCutValue/");
622 
623  cout<<"\n--------------------------------------------------------------------------------------" << endl;
624  Info("saveFWCutValuesHist()","Histogramms with CutValues for FW Hits saving:");
625  Int_t n=0;
626  for (int ep = 0; ep < (int) kNumFWCutValues; ++ep){
627  Int_t m = 0;
628  printf("\n>>>> %27s : ",getStringFWCutValues(ep).Data() );
629  for (int iCell = 0; iCell < (int) MAXFWCELLS; ++iCell){
630  if(fFWCutValuesHist[ep][iCell]){
631  TString temp = Form("FWCuts_%s_cell%d",getStringFWCutValues(ep).Data(), iCell );
632  fFWCutValuesHist[ep][iCell]->SetNameTitle(temp.Data(),temp.Data());
633  fFWCutValuesHist[ep][iCell]->Write();
634  if(m==0){
635  printf("%4s(%2d bins) ",fFWCutValuesHist[ep][iCell]->GetXaxis()->GetTitle(),
636  fFWCutValuesHist[ep][iCell]->GetNbinsX());
637  printf("%4s(%2d bins) ",fFWCutValuesHist[ep][iCell]->GetYaxis()->GetTitle(),
638  fFWCutValuesHist[ep][iCell]->GetNbinsY());
639  }
640  n++; m++;
641  }
642  }
643  printf("active cells:%2d ", m);
644 
645  }
646 
647  cout<< "\n\ntotal Histogramms with CutValues for FW Hits #"<< n << "saved" << endl;
648  cout<< "--------------------------------------------------------------------------------------" << endl;
649 
650  return 0;
651 }
652 
653 
655  return 0;
656 }
657 
659  return 0;
660 }
661 
663 {
664  // Reset.
665 }
666 //----------------------------------------------------------------------
668 {
670  return kFALSE;
671  }
672  else{
674  loadDayOfYear();
675  return kTRUE;
676  }
677 }
678 
680 {
681  if(isSimulation) return 1;
682 // ---- PT3 Event without PT2 Events --------------
683  if(gHades->getCurrentEvent()->getHeader()->getTBit() == 8192) return 1;
684 // ---- PT2 Event with PT3 Events --------------
685  if(gHades->getCurrentEvent()->getHeader()->getTBit() == 12288)return 1;
686 // ---- PT2 Event without PT3 Events ----------
687  if(gHades->getCurrentEvent()->getHeader()->getTBit() == 4096){
688  if(fDayOfYear<=106) return 4; //see logbook entry 15.4.2012 ~23.52h from MultSectra: 3.7;
689  else return 8; // 7.8;
690  }
691  return 0;
692 }
693 
694 
695 //----------------------------------------------------------------------
696 
698 {
699  Int_t i = wall->getCell();
700 
701  if ( i < 0 || i > 301 ) return -1;
702  else if ( i==65 || i==66 || i==77 || i==78 ) return -1; //beam hole
703  else if ( i<144 ) return 1; // small modules
704  else if ( i<208 ) return 2; // middle modules
705  else if ( i<302 ) return 3; // large modules
706  return -1;
707 }
708 
710 {
711  for (int ep = 0; ep < (int) kNumFWCutValues; ++ep){
712  if(useFWCut[ep] && !PassesCutsFW(wall, ep) ) return kFALSE;
713  }
714  return kTRUE;
715 }
716 
718 {
719  if(!wall) return kFALSE;
720  Int_t cell = wall->getCell();
721  if(cell<0 || cell>= MAXFWCELLS) return kFALSE;
722 
724  if(cell ==20 || cell ==53 || cell ==54
725  || cell ==64 || cell ==67 || cell ==79) return kFALSE; //noisy cells
726  }
727 
728  Float_t value =0;
729  if(eFWCut==kBetaCuts ) value = wall->getDistance() / wall->getTime()/ 299.792458;
730  if(eFWCut==kTimeCuts ) value = wall->getTime();
731  if(eFWCut==kChargeCuts) value = wall->getCharge();
732 
733 
734  if(fFWCutValuesHist[eFWCut][cell]){
735  if(fFWCutValuesHist[eFWCut][cell]->GetBinContent(fFWCutValuesHist[eFWCut][cell]->FindBin(value)) > 0) return kTRUE;
736  }
737  return kFALSE;
738 }
739 
741 {
742  Float_t Charge = 0;
743  Int_t Z = 0;
744 
745  if(!wall) return 0;
746  if(!isSimulation) Charge = wall->getCharge();
747  else Charge = 93.*pow(wall->getCharge(),0.46-0.006*sqrt(wall->getCharge())); // parametrization from R.Holzmann
748  Int_t cell = wall->getCell();
749  if(cell<0 || cell>= MAXFWCELLS) return 0;
750 
751  if(fFWCutValuesHist[kChargeCuts][cell]){
752  Z = fFWCutValuesHist[kChargeCuts][cell]->GetBinContent(fFWCutValuesHist[kChargeCuts][cell]->FindBin(Charge) );
753  }
754  if(Z>1) Z=Z-1;
755 
756  return Z;
757 }
758 
760 {
761  if(flag<kNumEventPlaneParameter) return (TH1F*)fEventPlaneCorrectionHist[flag];
762  else return 0;
763 }
764 
765 
767 {
768  Float_t fCorrection = 0.;
769  Double_t cent = getCentralityPercentile();
770  if(getEventPlaneCorrectionHist(flag))
771  fCorrection = getEventPlaneCorrectionHist(flag)->GetBinContent(getEventPlaneCorrectionHist(flag)->FindBin(cent, fDayOfYear) );
772  return fCorrection;
773 }
774 
776 {
777  Float_t fCorrectionError = 1.;
778  Double_t cent = getCentralityPercentile();
779  if(getEventPlaneCorrectionHist(flag))
780  fCorrectionError = getEventPlaneCorrectionHist(flag)->GetBinError(getEventPlaneCorrectionHist(flag)->FindBin(cent, fDayOfYear) );
781  return fCorrectionError;
782 }
783 
785  if(size==1) return ( fRandom->Rndm(1) - 0.5 )*40.;
786  else if(size==2) return ( fRandom->Rndm(1) - 0.5 )*80.;
787  else if(size==3) return ( fRandom->Rndm(1) - 0.5 )*160.;
788  return 0;
789 }
790 
791 Float_t HParticleEvtChara::getThetaWeight(HWallHitSim* wall, Float_t min, Float_t max){
792  Float_t theta = wall->getTheta();
793  if(theta>max) return 0.; // max Theta
794  if(theta<min) return 1.; // min Theta
795  Float_t wtheta = 1. -( (theta-min)/(max-min) );
796  if(wtheta>0 && wtheta < 1.) return wtheta;
797  else return 0.;
798 }
799 
801 {
802  // -----------------------------------------------------------------
803  //fill Hits into HitArray for QVector
804  fQVectorCalcDone = kFALSE; //reset QVectors
806  Int_t nHits = arrayOfHits.size();
807  for (Int_t i=0; i<nHits; i++){
808  if (!arrayOfHits[i]) continue;
809  delete arrayOfHits[i];
810  arrayOfHits[i] = NULL;
811  }
812  iFWHitvector.clear();
813  arrayOfHits.clear();
814 
815  if(!fCatWallHit) return;
816 
817  HWallHitSim *wall = 0; // includes HWallHit
818  Float_t wallXOrg=0,wallYOrg=0,wallZOrg=0;
819  Float_t wallX=0,wallY=0;
820  Float_t weight1 = 1.;
821  Float_t weight2 = 1.;
822  for(Int_t i = 0; i < fCatWallHit->getEntries(); i ++ ){
824  if(!PassesCutsFW(wall)) continue;
825  iFWHitvector.push_back(i);
826 
827  wall->getXYZLab(wallXOrg,wallYOrg,wallZOrg);
828 
829  wallX = (wallXOrg - getCorrection(kFWx)) / getCorrectionError(kFWx);
830  wallY = (wallYOrg - getCorrection(kFWy)) / getCorrectionError(kFWy);
831  weight1 = getFWCharge(wall);
832  weight2 = getThetaWeight(wall);
833 
834  if(wallX==0 || wallY==0) continue;
835  SimpleQVector* HitVector = new SimpleQVector();
836  HitVector->Set(wallX,wallY,weight1, weight2);
837  HitVector->SetOrg(wallXOrg,wallYOrg);
838  arrayOfHits.push_back(HitVector);
839 
840  }
841  nHits = arrayOfHits.size();
842  if(nHits<1) return;
843  //shuffle
844  std::random_shuffle ( arrayOfHits.begin(), arrayOfHits.end() );
845  //Flag two subevents in given shuffled sub-events
846  for (Int_t i=0; i<nHits; i++){
847  if (!arrayOfHits[i]) continue;
848  if((i+1)*2 > nHits) arrayOfHits[i]->SetSubEvent(1);
849  else arrayOfHits[i]->SetSubEvent(2);
850  }
851  return;
852 }
853 
855 {
856  // -----------------------------------------------------------------
857  if(isNewEvent()){
858  fillHitArray();
859  }
860  //print Hits of HitArray
861  Int_t nHits = arrayOfHits.size();
862  printf( "printHitArray #Hits = %d\n", nHits );
863  if(nHits<1) return;
864  for (Int_t i=0; i<nHits; i++){
865  if (arrayOfHits[i]) arrayOfHits[i]->print();
866  }
867  return;
868 }
869 
871  if(isNewEvent()){
872  fillHitArray();
873  }
874  return iFWHitvector;
875 }
876 
877 Bool_t HParticleEvtChara::fillQVectors(UInt_t statusflag, UInt_t nHarmonic)
878 {
879  // -----------------------------------------------------------------
880  //HWallHit into QVector
882  fQVectorCalcDone = kFALSE;
883 
884  vQPhi[0] = -1;
885  vQPhi[1] = -1;
886  vQPhi[2] = -1;
887 
888  Double_t dQX = 0.,dQY = 0.;
889  Double_t dQXA = 0.,dQYA = 0.;
890  Double_t dQXB = 0.,dQYB = 0.;
891 
892  Double_t dQXShift = 0.,dQYShift = 0.;
893  Double_t dQXScale = 1.,dQYScale = 1.;
894 
895  Double_t phi = -1;
896  Double_t w = 1.;
897 
898  Double_t sumOfWeights=0;
899 
900  Int_t nHits = arrayOfHits.size();
901  if(nHits<4) return kFALSE; // at least two hits per sub-event for
902  for (Int_t i=0; i<nHits; i++){
903  w = 1.;
904  if(isFlagSet(kShiftFW, statusflag)) phi = arrayOfHits[i]->Phi();
905  else phi = arrayOfHits[i]->PhiOrg();
906 
907  if(isFlagSet(kWeightCharge, statusflag)) w = arrayOfHits[i]->Weight1();
908  if(isFlagSet(kWeightTheta, statusflag)) w *= arrayOfHits[i]->Weight2();
909 
910  dQX += w * TMath::Cos(nHarmonic*phi);
911  dQY += w * TMath::Sin(nHarmonic*phi);
912 
913  if(arrayOfHits[i]->SubEvent() == 1) {
914  dQXA += w * TMath::Cos(nHarmonic*phi);
915  dQYA += w * TMath::Sin(nHarmonic*phi);
916  }
917  else if(arrayOfHits[i]->SubEvent() == 2) {
918  dQXB += w * TMath::Cos(nHarmonic*phi);
919  dQYB += w * TMath::Sin(nHarmonic*phi);
920  }
921  sumOfWeights +=w;
922  }
923 
924  if(isFlagSet(kReCentering, statusflag)){
925  if(isFlagSet(kShiftFW, statusflag)){
926  if(isFlagSet(kWeightCharge, statusflag)){ dQXShift = getCorrection(kQx2WCharge); dQYShift = getCorrection(kQy2WCharge);}
927  else{ dQXShift = getCorrection(kQx2); dQYShift = getCorrection(kQy2);}
928  }
929  else{
930  if(isFlagSet(kWeightCharge, statusflag)){ dQXShift = getCorrection(kQxWCharge); dQYShift = getCorrection(kQyWCharge);}
931  else{ dQXShift = getCorrection(kQx); dQYShift = getCorrection(kQy);}
932  }
933 
934  dQX -= dQXShift;
935  dQY -= dQYShift;
936  dQXA -= dQXShift;
937  dQYA -= dQYShift;
938  dQXB -= dQXShift;
939  dQYB -= dQYShift;
940  if(isFlagSet(kScaling, statusflag)){
941  if(isFlagSet(kShiftFW, statusflag)){
942  if(isFlagSet(kWeightCharge, statusflag)){dQXScale = getCorrectionError(kQx2WCharge); dQYScale = getCorrectionError(kQy2WCharge);}
943  else{ dQXScale = getCorrectionError(kQx2); dQYScale = getCorrectionError(kQy2);}
944  }
945  else{
946  if(isFlagSet(kWeightCharge, statusflag)){dQXScale = getCorrectionError(kQxWCharge); dQYScale = getCorrectionError(kQyWCharge);}
947  else{ dQXScale = getCorrectionError(kQx); dQYScale = getCorrectionError(kQy);}
948  }
949 
950  dQX /= dQXScale;
951  dQY /= dQYScale;
952  dQXA /= dQXScale;
953  dQYA /= dQYScale;
954  dQXB /= dQXScale;
955  dQYB /= dQYScale;
956  }
957  }
958 
959  if(dQX==0 && dQX==0) return kFALSE;
960  if(dQXA==0 && dQXA==0) return kFALSE;
961  if(dQXB==0 && dQXB==0) return kFALSE;
962 
963  if(isFlagSet(kRotation, statusflag)){
964  Float_t corrPsi = getCorrectionPhi(getPhi(dQX ,dQY ));
965  vQPhi[0] = TVector2::Phi_0_2pi( getPhi(dQX ,dQY ) + corrPsi );
966  vQPhi[1] = TVector2::Phi_0_2pi( getPhi(dQXA,dQYA) + corrPsi );
967  vQPhi[2] = TVector2::Phi_0_2pi( getPhi(dQXB,dQYB) + corrPsi );
968  }
969  else{
970  vQPhi[0] = getPhi(dQX ,dQY );
971  vQPhi[1] = getPhi(dQXA,dQYA);
972  vQPhi[2] = getPhi(dQXB,dQYB);
973  }
974 
975  fEventPlaneCorrectionFlag = statusflag;
976  fQVectorCalcDone = kTRUE;
977 
978  return kTRUE;
979 }
980 
982 {
983  // -----------------------------------------------------------------
984  if(fQVectorCalcDone){ printf("QVector Calculation Done\n");}
985  else{
986  printf("QVector Calculation NOT Done !!!!!!! retry now\n");
987  if(!fillQVectors(fEventPlaneCorrectionFlag)) {printf("QVector Calculation failed !!!!!!!\n"); return ;}
988  }
989 
990  printf("##### vQPhi: %.3f \t vQPhiA: %.3f \t vQPhiB: %.3f \t\t Corr: %s #####\n",
992  return;
993 }
994 
995 
997 {
998  //used Formular in Footnoot1 of Phys.Rev. C58 (1998) 1671-1678
999  //also appendix of Phys.Rev. C56 (1997) 3254-3264
1000 
1001  // kFourierCn is equivalent to Bn
1002  // kFourierSn is equivalent to An
1003 
1004 
1005  Float_t corrPsi = phi;
1006  corrPsi = 2./1.*(getCorrection(kFourierC1)*TMath::Sin(1.*phi))
1007  + 2./2.*(getCorrection(kFourierC2)*TMath::Sin(2.*phi))
1008  + 2./3.*(getCorrection(kFourierC3)*TMath::Sin(3.*phi))
1009  + 2./4.*(getCorrection(kFourierC4)*TMath::Sin(4.*phi))
1010  + 2./5.*(getCorrection(kFourierC5)*TMath::Sin(5.*phi))
1011  + 2./6.*(getCorrection(kFourierC6)*TMath::Sin(6.*phi))
1012  + 2./7.*(getCorrection(kFourierC7)*TMath::Sin(7.*phi))
1013  + 2./8.*(getCorrection(kFourierC8)*TMath::Sin(8.*phi))
1014 
1015  - 2./1.*(getCorrection(kFourierS1)*TMath::Cos(1.*phi))
1016  - 2./2.*(getCorrection(kFourierS2)*TMath::Cos(2.*phi))
1017  - 2./3.*(getCorrection(kFourierS3)*TMath::Cos(3.*phi))
1018  - 2./4.*(getCorrection(kFourierS4)*TMath::Cos(4.*phi))
1019  - 2./5.*(getCorrection(kFourierS5)*TMath::Cos(5.*phi))
1020  - 2./6.*(getCorrection(kFourierS6)*TMath::Cos(6.*phi))
1021  - 2./7.*(getCorrection(kFourierS7)*TMath::Cos(7.*phi))
1022  - 2./8.*(getCorrection(kFourierS8)*TMath::Cos(8.*phi));
1023  return corrPsi;
1024 }
1025 
1026 Float_t HParticleEvtChara::getEventPlane(UInt_t statusflag, UInt_t SubEvent, UInt_t nHarmonic)
1027 {
1028  if(isNewEvent()){
1029  fillHitArray();
1030  }
1031  if(!fQVectorCalcDone){
1032  if(!fillQVectors(statusflag)) return -1;
1033  }
1034 
1035  if(fQVectorCalcDone && fEventPlaneCorrectionFlag == statusflag){
1036  if(SubEvent<3) return vQPhi[SubEvent];
1037  else return -1;
1038  }
1039  else{
1040  if(!fillQVectors(statusflag, nHarmonic)) return -1;
1041  if(fQVectorCalcDone && fEventPlaneCorrectionFlag == statusflag){
1042  if(SubEvent<3) return vQPhi[SubEvent];
1043  else return -1;
1044  }
1045  }
1046  return -1;
1047 
1048 }
1049 
1050 Float_t HParticleEvtChara::getEventPlaneWeight(UInt_t statusflag, UInt_t SubEvent, UInt_t nHarmonic)
1051 {
1052  return 1.;
1053 }
1054 
1055 Float_t HParticleEvtChara::getEventPlaneParameter(UInt_t e, Bool_t corr)
1056 {
1057 
1058  if(isNewEvent()){
1059  fillHitArray();
1060  }
1061  if(e==kFWx|| e==kFWy) return -999;
1062 
1063  if(corr){
1064  if(e==kQx || e==kQy) fillQVectors(kReCentering|kScaling);
1069  }
1070  else{
1071  if(e==kQx || e==kQy) fillQVectors(kNoCorrection);
1073  if(e==kQx2 || e==kQy2|| e==kResolution) fillQVectors(kShiftFW);
1076  }
1077  if(vQPhi[0]==-1) return -999;
1078 
1079  if(e==kQx || e==kQxWCharge || e==kQx2 || e==kQx2WCharge) return TMath::Cos( vQPhi[0] );
1080  if(e==kQy || e==kQyWCharge || e==kQy2 || e==kQy2WCharge) return TMath::Sin( vQPhi[0] );
1081  if(e==kFourierC1) return TMath::Cos( vQPhi[0] );
1082  if(e==kFourierC2) return TMath::Cos(2.* vQPhi[0] );
1083  if(e==kFourierC3) return TMath::Cos(3.* vQPhi[0] );
1084  if(e==kFourierC4) return TMath::Cos(4.* vQPhi[0] );
1085  if(e==kFourierC5) return TMath::Cos(5.* vQPhi[0] );
1086  if(e==kFourierC6) return TMath::Cos(6.* vQPhi[0] );
1087  if(e==kFourierC7) return TMath::Cos(7.* vQPhi[0] );
1088  if(e==kFourierC8) return TMath::Cos(8.* vQPhi[0] );
1089  if(e==kFourierS1) return TMath::Sin( vQPhi[0] );
1090  if(e==kFourierS2) return TMath::Sin(2.* vQPhi[0] );
1091  if(e==kFourierS3) return TMath::Sin(3.* vQPhi[0] );
1092  if(e==kFourierS4) return TMath::Sin(4.* vQPhi[0] );
1093  if(e==kFourierS5) return TMath::Sin(5.* vQPhi[0] );
1094  if(e==kFourierS6) return TMath::Sin(6.* vQPhi[0] );
1095  if(e==kFourierS7) return TMath::Sin(7.* vQPhi[0] );
1096  if(e==kFourierS8) return TMath::Sin(8.* vQPhi[0] );
1097  if(e==kResolution || e==kResolutionWCharge){ //FIXME only 1 higher harmonics not included!!!
1098  return TMath::Cos(vQPhi[1])*TMath::Cos(vQPhi[2]) + TMath::Sin(vQPhi[1])*TMath::Sin(vQPhi[2]);
1099  }
1100  return -999;
1101 }
1102 
1103 
1104 //----------------------------------------------------------------------
1106 {
1107  // legacy code
1108  // Return centrality class, default in 5% of total cross section with estimator
1109  // or with preset classes like FOPI {6.3%, 21.0%, 30.9%}
1110 
1111  if (estimator.CompareTo("TOFRPCtimecut")==0 || estimator.CompareTo("TOFRPC5")==0)
1113  else if (estimator.CompareTo("TOFRPCtimecutFOPI")==0 || estimator.CompareTo("TOFRPCFOPI")==0)
1115  else if (estimator.CompareTo("TOFRPCtimecut10")==0 || estimator.CompareTo("TOFRPC10")==0)
1117  else {
1118  Error("getCentralityClass()","No CentralityEstimator defined!");
1119  return 0;
1120  }
1121 }
1122 
1124 {
1125  HParticleEvtInfo *event_info = (HParticleEvtInfo*)fParticleEvtInfoCat->getObject(0);
1126  if(!event_info) {
1127  Error("getCentralityClass()","No HParticleEvtInfo");
1128  return 0;
1129  }
1130 
1131  if(centE ==kTOFRPC) return event_info->getSumTofMultCut() + event_info->getSumRpcMultHitCut();
1132  else if(centE==kTOF) return event_info->getSumTofMultCut();
1133  else if(centE==kRPC) return event_info->getSumRpcMultHitCut();
1134  else if(centE==kTOFRPCtot) return event_info->getSumTofMult() + event_info->getSumRpcMultHit();
1135  else if(centE==kTOFtot) return event_info->getSumTofMult();
1136  else if(centE==kRPCtot) return event_info->getSumRpcMultHit();
1137  else if(centE==kSelectedParticleCand) return event_info->getSumSelectedParticleCandMult();
1138  else if(centE==kSelectedParticleCandCorr) return event_info->getSumSelectedParticleCandMult()*(fReferenceMeanSelTrack/event_info->getMeanMult());
1139  else if(centE==kSelectedParticleCandNorm) return event_info->getSumSelectedParticleCandMult()/event_info->getMeanMult();
1142  else if(centE==kPrimaryParticleCand) return event_info->getSumPrimaryParticleCandMult();
1143  else if(centE==kMdcWires) return event_info->getMdcWires();
1144  else if(centE==kMdcWiresOuterMod) return getMdcWiresOuterMod();
1145  else if(centE==kFWSumChargeSpec) return getFWSumChargeSpec();
1146  else if(centE==kFWSumChargeZ) return getFWSumZ();
1147  else if(centE==kDirectivity) return getDirectivity();
1148  else if(centE==kRatioEtEz) return getRatioEtEz();
1149  else if(centE==kEt) return getEt();
1150  else return 0;
1151 }
1152 
1153 TH1F* HParticleEvtChara::getCentralityClassHist(UInt_t centE, UInt_t centC) const
1154 {
1155  if(centE<kNumCentralityEstimator && centC<kNumCentralityClass) return (TH1F*)fCentralityHist[centE][centC];
1156  else return 0;
1157 }
1158 
1160 {
1161  if(centE<kNumCentralityEstimator) return (TH1F*)fCentralityPercentileHist[centE];
1162  else return 0;
1163 }
1164 
1165 Int_t HParticleEvtChara::getCentralityClass(UInt_t centE, UInt_t centC)
1166 {
1167  // Return centrality class, default in 5% of total cross section with estimator
1168  // or with preset classes like FOPI {6.3%, 21.0%, 30.9%}
1169  if(centE>=kNumCentralityEstimator){
1170  return 101.;
1171  }
1172  else if(!getCentralityClassHist(centE, centC)){
1173  //Error("getCentralityClass()","No CentralityEstimator defined!");
1174  return 0;
1175  }
1176  else {
1177  Int_t buffer = getCentralityEstimator(centE);
1178  Int_t bin = getCentralityClassHist(centE, centC)->FindBin(buffer);
1179  return getCentralityClassHist(centE, centC)->GetBinContent(bin);
1180  }
1181  return 0;
1182 
1183 }
1184 
1186 {
1187  // legacy code
1188  // print all CentralityClasses in the Estimator
1189 
1190  if (estimator.CompareTo("TOFRPCtimecut")==0 || estimator.CompareTo("TOFRPC5")==0)
1192  else if (estimator.CompareTo("TOFRPCtimecutFOPI")==0 || estimator.CompareTo("TOFRPCFOPI")==0)
1194  else if (estimator.CompareTo("TOFRPCtimecut10")==0 || estimator.CompareTo("TOFRPC10")==0)
1196  else {
1197  Error("printCentralityClass()","Sorry. printCentralityClass() for %s not implemented yet!",estimator.Data());
1198  return kFALSE;
1199  }
1200 
1201  return kFALSE;
1202 }
1203 
1204 Bool_t HParticleEvtChara::printCentralityClass(UInt_t centE, UInt_t centC)
1205 {
1206  // print all CentralityClasses in the Estimator
1207  cout << endl;
1208  cout<<"---------------------------------------------------------------------------------------------" << endl;
1209  cout<<"Centrality Classes for "<< getStringCentralityEstimator(centE).Data() << " with " << getStringCentralityClass(centC).Data() <<" bins" << endl;
1210  return printCentralityClass(getCentralityClassHist(centE, centC));
1211 }
1212 
1214 {
1215  // print all CentralityClasses in the Estimator
1216  if(!htemp){
1217  Error("printCentralityClass()","Sorry. printCentralityClass() for ... not implemented yet!");
1218  return kFALSE;
1219  }
1220  else{
1221  cout<<"---------------------------------------------------------------------------------------------" << endl;
1222  cout<< "# of Classes: "<< htemp->GetNbinsX()-2 << endl;
1223  printf(" Class lowEdge - upEdge Centrality[%%] BinWidth[%%] real CentralityBin[%%] BinCenter[%%]\n");
1224  Float_t pcent = 0;
1225  for(Int_t i = htemp->GetNbinsX(); i>0 ; --i) {
1226  printf(" %2.f : %8.2f - %8.2f %13s %13.3f %8.2f - %8.2f %13.2f\n",
1227  htemp->GetBinContent(i),
1228  (htemp->GetXaxis())->GetBinLowEdge(i),
1229  (htemp->GetXaxis())->GetBinUpEdge(i),
1230  (htemp->GetXaxis())->GetBinLabel(i),
1231  htemp->GetBinError(i),
1232  pcent,
1233  pcent+htemp->GetBinError(i),
1234  pcent+0.5*htemp->GetBinError(i) );
1235  pcent += htemp->GetBinError(i);
1236  }
1237  cout<<"---------------------------------------------------------------------------------------------" << endl;
1238  return kTRUE;
1239  }
1240  return kFALSE;
1241 }
1242 
1243 TH2D* HParticleEvtChara::makeEPresolution(TProfile2D *hist, Bool_t calcChi)
1244 {
1245  if(!hist) return 0;
1246  TString sName=hist->GetName();
1247  sName.Append("_res");
1248 
1249  Int_t nxBins = hist->GetNbinsX();
1250  Int_t nyBins = hist->GetNbinsY();
1251 
1252  TH2D* hresolution = (TH2D*) new TH2D(sName.Data(), sName.Data(), nxBins, hist->GetXaxis()->GetXmin(), hist->GetXaxis()->GetXmax(),
1253  nyBins, hist->GetYaxis()->GetXmin(), hist->GetYaxis()->GetXmax());
1254  (hresolution->GetXaxis())->SetTitle(hist->GetXaxis()->GetTitle());
1255  (hresolution->GetYaxis())->SetTitle(hist->GetYaxis()->GetTitle());
1256  (hresolution->GetZaxis())->SetTitle(hist->GetZaxis()->GetTitle());
1257 
1258  Int_t bin,binx,biny;
1259  for (binx =0;binx<=nxBins+1;binx++) {
1260  for (biny =0;biny<=nyBins+1;biny++) {
1261  bin = biny*(nxBins+2) + binx;
1262 
1263  Double_t dQaQb = hist->GetBinContent(bin);
1264  Double_t dEntriesQaQb = hist->GetBinEntries(bin);
1265  if( dQaQb <= 0 || dEntriesQaQb < 1 ){
1266  hresolution->SetBinContent(bin, 0. );
1267  hresolution->SetBinError( bin, 0. );
1268  continue;
1269  }
1270  Double_t dSpreadQaQb = hist->GetBinError(bin);
1271  Double_t dV = TMath::Sqrt(dQaQb);
1272  printf("\nbin=%d binX=%d binY=%d QaQb = %f +- %f V = %f\n",bin, binx, biny, dQaQb, dSpreadQaQb , dV);
1273  Double_t dChi = FindXi(dV,1e-6);
1274  dV = ComputeResolution( TMath::Sqrt2()*dChi );
1275  printf("An estimate of the event plane resolution is: %f\n", dV );
1276  Double_t dVerr = 0.;
1277  if(dQaQb > 0.) dVerr = (1./(2.*pow(dQaQb,0.5)))*dSpreadQaQb;
1278  Double_t dChiErr = FindXi(dVerr,1e-6);
1279  printf("An estimate chi of the event plane is: %f +- %f\n", dChi, dChiErr );
1280  printf("R:(subevents) = %f +- %f\n",dV,dVerr);
1281  if(calcChi){
1282  hresolution->SetBinContent(binx, biny, dChi );
1283  hresolution->SetBinError( binx, biny, dChiErr);
1284  }
1285  else{
1286  hresolution->SetBinContent(binx, biny, dV );
1287  hresolution->SetBinError( binx, biny, dVerr );
1288  }
1289  }
1290  }
1291  TString name=hist->GetName();
1292  hresolution->SetNameTitle(name.Data(),name.Data());
1293  hresolution->SetEntries(hist->GetEntries());
1294  return hresolution;
1295 
1296 }
1297 TH1* HParticleEvtChara::makeEPresolution(TH3 *hist, Bool_t calcChi)
1298 {
1299  if(!hist) return 0;
1300  TString sName=hist->GetName();
1301  sName.Append("_res");
1302 
1303  TAxis* axis = hist->GetZaxis();
1304  axis->SetRange(1, axis->GetNbins());
1305  TH1* temp = hist->Project3D("xy_1");
1306 
1307  axis->SetRange(axis->FindBin(TMath::PiOver2()), axis->GetNbins());
1308  TH1* temp2 = hist->Project3D("xy_2");
1309  temp2->Divide(temp);
1310 
1311  Int_t nxBins = temp2->GetNbinsX();
1312  Int_t nyBins = temp2->GetNbinsY();
1313 
1314  TH2D* hresolution = (TH2D*) new TH2D(sName.Data(), sName.Data(), nxBins, temp2->GetXaxis()->GetXmin(), temp2->GetXaxis()->GetXmax(),
1315  nyBins, temp2->GetYaxis()->GetXmin(), temp2->GetYaxis()->GetXmax());
1316  (hresolution->GetXaxis())->SetTitle(hist->GetXaxis()->GetTitle());
1317  (hresolution->GetYaxis())->SetTitle(hist->GetYaxis()->GetTitle());
1318  (hresolution->GetZaxis())->SetTitle(hist->GetZaxis()->GetTitle());
1319 
1320  Int_t bin,binx,biny;
1321 
1322  for (binx =0;binx<=nxBins+1;binx++) {
1323  for (biny =0;biny<=nyBins+1;biny++) {
1324  bin = biny*(nxBins+2) + binx;
1325 
1326  Double_t ratio = temp2->GetBinContent(binx,biny);
1327  printf("\nbin=%d binX=%d binY=%d Ratio=%f\n",bin, binx, biny, ratio);
1328  if( ratio <= 0 ){
1329  continue;
1330  }
1331  Double_t chisq = -2.*TMath::Log(2.*ratio);
1332  Double_t dChi = sqrt(chisq);
1333  printf("An estimate chi of the event plane is: %f\n", dChi );
1334  Double_t dV = ComputeResolution( TMath::Sqrt2()*dChi );
1335  printf("An estimate of the event plane resolution is: %f\n", dV );
1336  Double_t dVerr = 0.;
1337  printf("R:(subevents) = %f +- %f\n",dV,dVerr);
1338 
1339  if(calcChi){
1340  hresolution->SetBinContent(binx, biny, dChi );
1341  //hresolution->SetBinError( binx, biny, dChiErr);
1342  }
1343  else{
1344 
1345  hresolution->SetBinContent(binx, biny, dV );
1346  hresolution->SetBinError( binx, biny, dVerr );
1347  }
1348  }
1349  }
1350 
1351  TString name=hist->GetName();
1352  hresolution->SetNameTitle(name.Data(),name.Data());
1353  hresolution->SetEntries(hist->GetEntries());
1354  return hresolution;
1355 }
1356 
1357 //--------------------------------------------------------------------
1358 Double_t HParticleEvtChara::ModifiedBesselI(Int_t n, Double_t x) const
1359 {
1360  // compute half-integer modified Bessel functions
1361  // order: n>0, for n>5/2, interpolation is used (improve this by using recurrence!!!)
1362  const Double_t FACTOR = 0.797884561; //FIXME
1363  if (n<0) return 0;
1364  if (x<1e-7) return 0;
1365 
1366  if (n==0) return FACTOR*sqrt(x)*TMath::SinH(x)/x; // 1/2
1367  else if (n==1) return FACTOR*sqrt(x)*( -TMath::SinH(x)/(x*x) + TMath::CosH(x)/x ); // 3/2
1368  else if (n==2) return FACTOR*sqrt(x)*( (3./(x*x)+1.)*TMath::SinH(x)/x - 3.*TMath::CosH(x)/(x*x) ); // 5/2
1369  return 0.5*(TMath::BesselI(n,x)+TMath::BesselI(n+1,x)); // use average of integer-order Bessel
1370 }
1371 
1372 //--------------------------------------------------------------------
1373 Double_t HParticleEvtChara::ComputeResolution( Double_t x, Int_t n ) const
1374 {
1375  // Computes resolution for Event Plane method
1376  if(x > 51.3) {
1377  printf("Warning: Estimation of total resolution might be WRONG. Please check!");
1378  return 0.99981;
1379  }
1380  if (n<1) return 0;
1381  Int_t n1= (n-1)/2;
1382  Int_t n2 = (n+1)/2;
1383 
1384  Double_t a = x*x/2; // in formula (8) of Ollitrault arXiv:nucl-ex/9711003v2
1385  // it is x*x/2 for full EP resolution here sub-event resolution
1386  Double_t b = TMath::Exp(-a);
1387  if (n==1) b *= TMath::BesselI0(a)+TMath::BesselI1(a);
1388  else if(n%2==1) b *= TMath::BesselI(n1,a)+TMath::BesselI(n2,a);
1389  else b *= ModifiedBesselI(n1, a) + ModifiedBesselI(n2, a);
1390  return TMath::Sqrt(TMath::Pi())/2*x*b;
1391 }
1392 
1393 
1394 //--------------------------------------------------------------------
1395 Double_t HParticleEvtChara::FindXi( Double_t res, Double_t prec, Int_t n ) const
1396 {
1397  // Computes x(res) for Event Plane method
1398  if(res > 0.99981) {
1399  printf("Warning: Resolution for subEvent is high. You reached the precision limit.");
1400  return 51.3;
1401  }
1402  int nSteps =0;
1403  Double_t xtmp=0, xmin=0, xmax=51.3, rtmp=0, delta=2*prec;
1404  while( delta > prec ) {
1405  xtmp = 0.5*(xmin+xmax);
1406  rtmp = ComputeResolution(xtmp, n);
1407  delta = TMath::Abs( res-rtmp );
1408  if(rtmp>res) xmax = xtmp;
1409  if(rtmp<res) xmin = xtmp;
1410  nSteps++;
1411  }
1412  return xtmp;
1413 }
1414 //--------------------------------------------------------------------
1415 Bool_t HParticleEvtChara::addEstimatorHist(TH1F *hist, Float_t fractionXsection, UInt_t centE, Int_t direction)
1416 {
1417  if(!hist)return kFALSE;
1418  if(fractionXsection<0.1) return kFALSE;
1419  if(centE>=kNumCentralityEstimator) return kFALSE;
1420 
1421  TString temp;
1422  if(fPostFix.CompareTo("")==0) temp = Form("%s",getStringCentralityEstimator(centE).Data() );
1423  else temp = Form("%s_%s",getStringCentralityEstimator(centE).Data(), fPostFix.Data() );
1424 
1425  hist->SetName(temp.Data());
1426  hist->SetTitle(temp.Data());
1427  fEstimatorHist[centE]= (TH1F*) hist;
1428  fCentralityPercentileHist[centE] = (TH1F*) makePercentiles(hist, fractionXsection, direction);
1429  for (int centC = 0; centC < (int) kNumCentralityClass; ++centC){
1430  fCentralityHist[centE][centC] = (TH1F*) makeClasses(hist, fractionXsection, centC, direction);
1431  if(fCentralityHist[centE][centC]) printCentralityClass(centE, centC);
1432  }
1433  return kTRUE;
1434 }
1435 
1436 TH1F* HParticleEvtChara::makePercentiles(TH1F *htemp, Float_t fractionXsection, Int_t direction)
1437 {
1438  if(!htemp) return 0;
1439 
1440  TH1F *hpercent = (TH1F*) htemp->Clone("hpercent");
1441  TString name=htemp->GetName();
1442  name.Append("_percentile");
1443  hpercent->SetNameTitle(name.Data(),name.Data());
1444  hpercent->Reset();
1445  Float_t totIntegral = fractionXsection / htemp->Integral(1,htemp->GetNbinsX());
1446  if(direction<0){
1447  for (int ibin=1; ibin<=htemp->GetNbinsX(); ibin++){
1448  hpercent->SetBinContent(ibin, totIntegral * htemp->Integral(ibin,htemp->GetNbinsX()) );
1449  }
1450  }
1451  else if(direction>0){
1452  for (int ibin=1; ibin<=htemp->GetNbinsX(); ibin++){
1453  hpercent->SetBinContent(ibin, totIntegral * htemp->Integral(1,ibin) );
1454  }
1455  }
1456  return hpercent;
1457 }
1458 
1459 TH1F* HParticleEvtChara::makeClasses(TH1F *htemp, Float_t fractionXsection, UInt_t centC, Int_t direction)
1460 {
1461 
1462  if(!htemp) return 0;
1463  if(fractionXsection<0.1) return 0;
1464  Int_t nClasses = getCentralityClassNbins(centC);
1465  Double_t integral = htemp->Integral();
1466  Double_t norm = integral / fractionXsection;
1467  Float_t* PercentileArray = getCentralityClassArray(centC);
1468  std::vector<Double_t> binEdge;
1469  std::vector<Double_t> xSection;
1470  std::vector<TString> fLabels;
1471 
1472 
1473  Int_t start = 1;
1474  Int_t stop = htemp->GetNbinsX();
1475  if(direction < 0){
1476  direction = -1;
1477  start = htemp->GetNbinsX()-1;
1478  stop = 0;
1479  binEdge.push_back(htemp->GetBinLowEdge(htemp->GetNbinsX()));
1480  }
1481  else{
1482  direction = 1;
1483  binEdge.push_back(htemp->GetBinLowEdge(1));
1484  }
1485  Int_t bin = start;
1486  Float_t lxs = 0;
1487  Float_t txs = 0;
1488 
1489 
1490  // Find edge and starting point at begining
1491  while(1){
1492  lxs += htemp->GetBinContent(bin);
1493  Double_t pxs = lxs/norm;
1494  if( pxs>0.001 ){
1495  if (direction>0) binEdge.push_back(htemp->GetBinLowEdge(bin+1));
1496  else binEdge.push_back(htemp->GetBinLowEdge(bin));
1497  xSection.push_back(pxs);
1498  bin += direction;
1499  break;
1500  }
1501  bin += direction;
1502  }
1503 
1504  Float_t totxs = 0;
1505  for (Int_t i = 0; i < nClasses; ++i) {
1506  lxs = 0 ;
1507  totxs += PercentileArray[i];
1508  while(1){
1509  lxs += htemp->GetBinContent(bin);
1510  txs += htemp->GetBinContent(bin);
1511  Double_t pxs = txs/norm;
1512  Double_t tdiff = (txs+htemp->GetBinContent(bin+direction))/norm;
1513  if( ( pxs>totxs )
1514  || (TMath::Abs(pxs-totxs)<=TMath::Abs(tdiff-totxs)) )
1515  {
1516  if (direction>0) binEdge.push_back(htemp->GetBinLowEdge(bin+1));
1517  else binEdge.push_back(htemp->GetBinLowEdge(bin));
1518  xSection.push_back(lxs/norm);
1519  bin += direction;
1520  break;
1521  }
1522  bin += direction;
1523  if ( (direction>0 && bin>=stop) || (direction<0 && bin<=stop ) || (txs>=integral) ) break;
1524  }
1525  if(totxs>fractionXsection || (direction>0 && bin>=stop) || (direction<0 && bin<=stop ) || (txs>=integral)) break;
1526  }
1527 
1528  fLabels.push_back("overflow");
1529  Double_t totXsection =0;
1530  for(std::vector<double>::size_type index = 1; index < xSection.size()-1; ++index){
1531  fLabels.push_back(Form("%02.0f-%02.0f",round(totXsection),round(totXsection+xSection[index])) );
1532  totXsection += xSection[index];
1533  }
1534  fLabels.push_back("underflow");
1535 
1536  std::reverse(fLabels.begin(),fLabels.end());
1537  std::reverse(binEdge.begin(), binEdge.end());
1538  std::reverse(xSection.begin(), xSection.end());
1539 
1540  Double_t xlowbins[binEdge.size()];
1541  std::copy(binEdge.begin(), binEdge.end(), xlowbins);
1542 
1543  TString name = htemp->GetTitle();
1544  name = Form("%s_%s_fixedCuts", name.Data(), getStringCentralityClass(centC).Data());
1545  TH1F *hfixedCuts = new TH1F(name.Data(), name.Data(), binEdge.size()-1, xlowbins);
1546 
1547  for(std::vector<double>::size_type bin = 0; bin < fLabels.size(); ++bin){
1548  (hfixedCuts->GetXaxis())->SetBinLabel(bin+1,fLabels[bin]);
1549  hfixedCuts->SetBinContent(bin+1, fLabels.size()-bin-1);
1550  hfixedCuts->SetBinError(bin+1, xSection[bin]);
1551  }
1552  return hfixedCuts;
1553 }
1554 
1555 
1556 
1557 Int_t HParticleEvtChara::getNbins(TString estimator)
1558 {
1559  // legacy code
1560  // Number of Bins (CentralityClasses and Over- and Underflow) in the Estimator
1561 
1562  if (estimator.CompareTo("TOFRPCtimecut")==0 || estimator.CompareTo("TOFRPC5")==0)
1564  else if (estimator.CompareTo("TOFRPCtimecutFOPI")==0 || estimator.CompareTo("TOFRPCFOPI")==0)
1566  else if (estimator.CompareTo("TOFRPCtimecut10")==0 || estimator.CompareTo("TOFRPC10")==0)
1568  else { Error("getNbins()","Sorry. getNbins() for %s not implemented yet!",estimator.Data()); return 0;}
1569 }
1570 Int_t HParticleEvtChara::getNbins(UInt_t centE, UInt_t centC)
1571 {
1572  // Number of Bins (CentralityClasses and Over- and Underflow) in the Estimator
1573  TH1F *htemp = getCentralityClassHist(centE, centC);
1574  if(htemp) return htemp->GetNbinsX();
1575  else return 0;
1576 }
1577 
1579 {
1580  if(centC==k1040) return 6;
1581  else if(centC==kFOPI) return 6;
1582  Float_t binSize = getCentralityClassBinSize(centC);
1583  if(binSize>0) return round((100./binSize));
1584  // cout << "binSize: "<< binSize <<" nBins: " <<nBins << " " <<endl;
1585  return 0;
1586 }
1587 
1589 {
1590  //
1591  Float_t binSize = getCentralityClassBinSize(centC);
1592  Int_t nBins = getCentralityClassNbins(centC);
1593  if (nBins<1) return 0;
1594 
1595  Float_t* arr = new Float_t[nBins];
1596  if(binSize>0){
1597  Float_t xs = 0.;
1598  for(Int_t i =0; i < nBins ; i++){
1599  xs+=binSize;
1600  if(xs < 100.){
1601  arr[i] = (Float_t) binSize;
1602  }
1603  else{
1604  arr[i] = 100. - (binSize*i);
1605  }
1606  }
1607  }
1608  else if(centC==k1040) {
1609  Double_t fxs[6] = {10.,30.,30, 10.,10.,10.};
1610  for(Int_t i =0; i < nBins ; i++) arr[i] = fxs[i];
1611  }
1612  else if(centC==kFOPI) {
1613  Double_t fxs[6] = {6.3,14.7,9.9,10.,10.,10.};
1614  for(Int_t i =0; i < nBins ; i++) arr[i] = fxs[i];
1615  }
1616 
1617  return arr;
1618 }
1619 
1620 Float_t* HParticleEvtChara::getUpEdgeArray(UInt_t centE, UInt_t centC)
1621 {
1622  //
1623  TH1F *htemp = getCentralityClassHist(centE, centC);
1624 
1625  if(htemp){
1626  Int_t nBins = htemp->GetNbinsX();
1627  Float_t* arr = new Float_t[nBins];
1628  for(Int_t i =0; i < nBins ; i++) {
1629  //printf(" %8.2f\n",(htemp->GetXaxis())->GetBinLowEdge(nBins-i));
1630  arr[i] = (Float_t) (htemp->GetXaxis())->GetBinLowEdge(nBins-i);
1631  }
1632  return arr;
1633  }
1634  else return 0;
1635 }
1636 
1637 Float_t* HParticleEvtChara::getBinCenterArray(UInt_t centE, UInt_t centC)
1638 {
1639  //
1640  TH1F *htemp = getCentralityClassHist(centE, centC);
1641 
1642  if(htemp){
1643  Int_t nBins = htemp->GetNbinsX();
1644  Float_t* arr = new Float_t[nBins];
1645  Float_t pcent = 0;
1646  for(Int_t i =0; i < nBins ; i++) {
1647  pcent += 0.5*htemp->GetBinError(i);
1648  arr[nBins-i] = (Float_t) pcent;
1649  }
1650  return arr;
1651  }
1652  else return 0;
1653 }
1654 
1655 vector<TString> HParticleEvtChara::getLabelArray(UInt_t centE, UInt_t centC)
1656 {
1657  //
1658  TH1F *htemp = getCentralityClassHist(centE, centC);
1659  vector<TString> ar;
1660  if(htemp){
1661  Int_t nBins = htemp->GetNbinsX();
1662  for(Int_t i =1; i < nBins-1 ; i++) {
1663  //printf(" %s\n",(htemp->GetXaxis())->GetBinLabel(nBins-i));
1664  TString prv = TString( (htemp->GetXaxis())->GetBinLabel(nBins-i) );
1665  ar.push_back(prv);
1666  }
1667  }
1668  return ar;
1669 }
1670 
1672 {
1673  // Return centrality class 5% of total cross section
1674  return 1+(Int_t)getCentralityPercentile(estimator)/5.;
1675 }
1677 {
1678  // Return centrality class 5% of total cross section
1679  return 1+(Int_t)getCentralityPercentile(estimator)/10.;
1680 }
1681 
1683 {
1684  // legacy code
1686 }
1687 
1688 
1690 {
1691 
1692  if(centE>=kNumCentralityEstimator){
1693  return 101.;
1694  }
1695  else if(!getCentralityPercentileHist(centE)){
1696  return 101.;
1697  }
1698  else {
1699  Int_t buffer = getCentralityEstimator(centE);
1700  Int_t bin = getCentralityPercentileHist(centE)->FindBin(buffer);
1701  return getCentralityPercentileHist(centE)->GetBinContent(bin);
1702  }
1703 
1704  return 101.;
1705 }
1706 
1708 {
1709  HParticleEvtInfo *event_info = (HParticleEvtInfo*)fParticleEvtInfoCat->getObject(0);
1710  if(!event_info) {
1711  return 0;
1712  }
1713  Int_t sum = 0;
1714  for(Int_t s=0;s<6;s++) sum += event_info->getMdcWiresMod(s,2)+event_info->getMdcWiresMod(s,3); //outer mdc
1715  return sum;
1716 }
1717 
1718 
1720 
1721  if(!gLoop) return 0;
1722 
1723  HParticleEvtInfo *event_info = (HParticleEvtInfo*)fParticleEvtInfoCat->getObject(0);
1724  if(!event_info) {
1725  return 0;
1726  }
1727 
1728  Float_t SelectedParticleCandSecCorr=0;
1729  Int_t nGoodSectors=0;
1730  for (Int_t s = 0; s < 6; ++s ){
1731  if(gLoop->goodSector(s)){
1732  SelectedParticleCandSecCorr += event_info->getSelectedParticleCandMult(s)/event_info->getMeanMult(s);
1733  nGoodSectors++;
1734  }
1735  }
1736  if(nGoodSectors>0) SelectedParticleCandSecCorr /= nGoodSectors;
1737  return SelectedParticleCandSecCorr;
1738 }
1739 
1741 
1742  HParticleEvtInfo *event_info = (HParticleEvtInfo*)fParticleEvtInfoCat->getObject(0);
1743  if(!event_info) {
1744  return 0;
1745  }
1746  // Correction function with parameters provided by O.&V.Pechenova
1747  Float_t w0 = -0.660682;
1748  Float_t w1 = -0.0652827;
1749  Float_t w2 = -0.660682;
1750  Float_t w3 = -0.0652827;
1751  Float_t w4 = -0.655548;
1752  Float_t w5 = -0.00547515;
1753 
1754  Int_t nRawWires = 0;
1755  Int_t nSelMult = 0;
1756  Float_t WirePerTracks = 0;
1757  Float_t SelectedParticleCandCorrPerWire = 0;
1758  Float_t EffWirePerTrack = 0;
1759  for(Int_t s=0;s<6;s++) {
1760  nSelMult = event_info->getSelectedParticleCandMult(s);
1761  nRawWires = event_info->getMdcWiresMod(s,2)+event_info->getMdcWiresMod(s,3); //outer mdc
1762  if(nSelMult>0){
1763  WirePerTracks = Float_t(nRawWires)/nSelMult;
1764  EffWirePerTrack = 1/(exp(w0+w1*WirePerTracks) + exp(w2+w3*WirePerTracks) + exp(w4+w5*WirePerTracks));
1765  SelectedParticleCandCorrPerWire += nSelMult*EffWirePerTrack;
1766  }
1767  }
1768  return SelectedParticleCandCorrPerWire;
1769 }
1770 
1772 
1773  Float_t MultFWSumChargeSpec=0;
1774 
1775  //HWallHit
1776  if(!fCatWallHit) return -1;
1777  HWallHitSim *wall = 0;
1778  for(Int_t i = 0; i < fCatWallHit->getEntries(); i ++ ){
1780  if(!wall) continue;
1781  if(!PassesCutsFW(wall)) continue;
1782  if(!isSimulation) MultFWSumChargeSpec += wall->getCharge();
1783  else MultFWSumChargeSpec += 93.*pow(wall->getCharge(),0.46-0.006*sqrt(wall->getCharge())); // parametrization from R.Holzmann
1784  }
1785 
1786  return MultFWSumChargeSpec;
1787 }
1788 
1789 Float_t HParticleEvtChara::getFWSumZ(Int_t minZ, Int_t maxZ, UInt_t SubEvent){
1790 
1791  if(isNewEvent()){
1792  fillHitArray();
1793  }
1794  Float_t Z = 0.;
1795  Float_t w = 1.;
1796 
1797  Int_t nHits = arrayOfHits.size();
1798  if(nHits<1) return 0; // at least one hit
1799  for (Int_t i=0; i<nHits; i++){
1800  w = arrayOfHits[i]->Weight1();
1801  //w2 *= arrayOfHits[i]->Weight2();
1802  if(w>maxZ) w=maxZ;
1803  if(w<minZ) continue;
1804 
1805  if(SubEvent == 0) {
1806  Z += w;
1807  }
1808  else if(SubEvent == 1 && arrayOfHits[i]->SubEvent() == 1) {
1809  Z += w;
1810  }
1811  else if(SubEvent == 2 && arrayOfHits[i]->SubEvent() == 2) {
1812  Z += w;
1813  }
1814  }
1815  return Z;
1816 }
1817 
1819 
1820  Float_t Et=0;
1821  if(!fParticleCandCat) return -1;
1822  HParticleCand* particle_cand = 0;
1823 
1824  for(Int_t i = 0; i < fParticleCandCat->getEntries(); i ++ ){
1825  particle_cand = (HParticleCand*) HCategoryManager::getObject(particle_cand,fParticleCandCat,i);
1826  if(!particle_cand->isFlagBit(Particle::kIsUsed)) continue;
1827  if(particle_cand->isFakeRejected() ) continue;
1828  if(particle_cand->getSystemUsed() == -1) continue;
1829  if(!particle_cand->isFlagAND(4,
1833  Particle::kIsAcceptedRK) ) continue;
1834  if(particle_cand->getPID()==-1) continue;
1835 // if( particle_cand->getMomentumOrg() == particle_cand->getMomentum() ) continue;
1836  if( TMath::Abs(particle_cand->Rapidity()-0.74)>0.5 && particle_cand->Pt()>300 ){ //
1837  Et += particle_cand->Et();
1838  }
1839  //-----------------------------------------------------------------
1840  } // end particle loop
1841  return Et;
1842 }
1843 
1845 
1846  Float_t RatioEtEz=0;
1847  Float_t Et=0;
1848  Float_t Ez=0;
1849  if(!fParticleCandCat) return -1;
1850  HParticleCand* particle_cand = 0;
1851 
1852  for(Int_t i = 0; i < fParticleCandCat->getEntries(); i ++ ){
1853  particle_cand = (HParticleCand*) HCategoryManager::getObject(particle_cand,fParticleCandCat,i);
1854  if(!particle_cand->isFlagBit(Particle::kIsUsed)) continue;
1855  if(particle_cand->isFakeRejected() ) continue;
1856  if(particle_cand->getSystemUsed() == -1) continue;
1857  if(!particle_cand->isFlagAND(4,
1861  Particle::kIsAcceptedRK) ) continue;
1862  if(particle_cand->getPID()==-1) continue;
1863  Ez += particle_cand->E()*particle_cand->CosTheta(); // pt2 == 0 ? 0 : E()*E() * pt2/(pt2+Z()*Z());
1864  Et += particle_cand->Et();
1865  //-----------------------------------------------------------------
1866  } // end particle loop
1867  if(Ez){
1868  RatioEtEz = Et/Ez;
1869  }
1870  else {
1871  RatioEtEz=0;
1872  }
1873  return RatioEtEz;
1874 }
1875 
1876 
1878 {
1879  Float_t Directivity=0;
1880  TLorentzVector QVector = TLorentzVector();
1881  Int_t QVectorT = 0;
1882  if(!fParticleCandCat) return -1;
1883  HParticleCand* particle_cand = 0;
1884 
1885  for(Int_t i = 0; i < fParticleCandCat->getEntries(); i ++ ){
1886  particle_cand = (HParticleCand*) HCategoryManager::getObject(particle_cand,fParticleCandCat,i);
1887  if(!particle_cand->isFlagBit(Particle::kIsUsed)) continue;
1888  if(particle_cand->isFakeRejected() ) continue;
1889  if(particle_cand->getSystemUsed() == -1) continue;
1890  if(!particle_cand->isFlagAND(4,
1894  Particle::kIsAcceptedRK) ) continue;
1895  if(particle_cand->getPID()==-1) continue;
1896  if( particle_cand->Rapidity() > 0.74 ){ //
1897  QVector += (TLorentzVector)* particle_cand;
1898  QVectorT += particle_cand->Pt();
1899  }
1900  else { //
1901  QVector -= (TLorentzVector)* particle_cand;
1902  QVectorT += particle_cand->Pt();
1903  }
1904  //-----------------------------------------------------------------
1905  } // end particle loop
1906  if(QVectorT>0)Directivity = QVector.Pt() / QVectorT; // = |Sum(vector_p_{t,i})| / Sum(|p_{t,i}|) for 0>Y_cm
1907  else Directivity = 0;
1908  return Directivity;
1909 }
1910 
TString getStringCentralityClass(UInt_t e)
Float_t getCorrectionError(UInt_t flag=kQx)
Float_t * getCentralityClassArray(UInt_t centC=k10)
R__EXTERN HLoop * gLoop
Definition: hloop.h:134
Int_t getFWCharge(HWallHitSim *wall)
Int_t GetFWmoduleSize(HWallHitSim *wall)
UInt_t getEventSeqNumber(void)
Definition: heventheader.h:132
vector< TString > getLabelArray(UInt_t centE=kTOFRPC, UInt_t centC=k10)
TString getStringCentralityEstimator(UInt_t e)
Double_t ComputeResolution(Double_t x, Int_t n=1) const
Float_t * getBinCenterArray(UInt_t centE=kTOFRPC, UInt_t centC=k10)
Float_t getSelectedParticleCandSecNorm()
Int_t getNbins(TString estimator)
vector< Int_t > getFWhits()
Float_t getDistance(void) const
Definition: hwallhit.h:41
Float_t getEventPlane(UInt_t EPcorr=kNoCorrection, UInt_t SubEvent=0, UInt_t nHarmonic=1)
Float_t getCentralityClassBinSize(UInt_t e)
Bool_t isFlagBit(eFlagBits bit)
Bool_t PassesCutsFW(HWallHitSim *wall)
Bool_t addEstimatorHist(TH1F *hist, Float_t fractionXsection=100., UInt_t centE=kTOFRPC, Int_t direction=-1)
Bool_t saveCentralityEstimatorHist()
Double_t FindXi(Double_t res, Double_t prec=1e-6, Int_t n=1) const
Float_t getSelectedParticleCandCorrPerWire()
UInt_t getCentralityEstimatorFromString(TString s)
HDataSource * getDataSource(void) const
Definition: hades.h:124
vector< TH1 * > fEstimatorHist
TH1F * getCentralityPercentileHist(UInt_t centE=kTOFRPC) const
Int_t getSumRpcMultHit() const
HRuntimeDb * getRuntimeDb(void)
Definition: hades.h:111
Bool_t saveEventPlaneCorrectionHist()
HEvent *& getCurrentEvent(void)
Definition: hades.cc:422
const Cat_t catGeantKine
Definition: hgeantdef.h:8
Float_t getVersion()
destructor
TString getStringEventPlaneCorrection(UInt_t e)
p1 SetTitle("Momentum")
vector< Bool_t > useFWCut
virtual HEventHeader * getHeader(void) const
Definition: hevent.h:36
Float_t getCorrection(UInt_t flag=kQx)
Int_t n
Float_t getThetaWeight(HWallHitSim *wall, Float_t min=3.5, Float_t max=8.)
void Set(Float_t x, Float_t y, Float_t w1=1., Float_t w2=1.)
Int_t getMdcWiresMod(Int_t s, Int_t m)
Double_t theta
p1 GetXaxis() -> SetLabelSize(0.06)
static T * getObject(T *pout, Short_t num=-1, Int_t index=-1, Bool_t silent=kFALSE)
Float_t getPhi(Float_t x, Float_t y)
static HCategory * getCategory(Short_t num, Int_t warn=0, TString name="")
Bool_t isFakeRejected(Int_t io=-1)
Text_t const * getCurrentFileName()
Definition: hruntimedb.h:53
Bool_t loadCentralityEstimatorHist()
HCategory * fParticleEvtInfoCat
Float_t getTheta(void) const
Definition: hwallhit.h:39
ClassImp(HParticleEvtChara) HParticleEvtChara
Float_t * getUpEdgeArray(UInt_t centE=kTOFRPC, UInt_t centC=k10)
p1 GetYaxis() -> SetLabelSize(0.06)
Bool_t addFWCutValuesHist(TH1 *hist, Int_t cell, UInt_t eFWCut)
TH1F * makePercentiles(TH1F *hist, Float_t fractionXsection=100., Int_t direction=-1)
Float_t getCentralityPercentile(TString estimator)
Double_t sum
Int_t getCell(void) const
Definition: hwallhit.h:37
TH1F * makeClasses(TH1F *h, Float_t fractionXsection=100., UInt_t centC=k10, Int_t direction=-1)
Int_t getCentralityClass(TString estimator)
Int_t getSumSelectedParticleCandMult() const
Bool_t setParameterFile(TString ParameterFile)
Float_t getCorrectionPhi(Float_t phi)
void SetOrg(Float_t x, Float_t y)
HCategory * fParticleCandCat
vector< TH2 * > fEventPlaneCorrectionHist
static Int_t getDayFileName(TString name, Bool_t print=kFALSE)
Definition: htime.cc:170
const Cat_t catWallHit
Definition: walldef.h:9
~HParticleEvtChara()
constructor
Float_t getCentralityEstimator(UInt_t centE=kTOFRPC)
Hades * gHades
Definition: hades.cc:1213
Int_t getCentralityClass5(TString estimator)
static TString stripFileName(TString name, Bool_t removeEvtBuilder=kFALSE, Bool_t silent=kFALSE)
Definition: htime.cc:54
TH1F * getCentralityClassHist(UInt_t centE=kTOFRPC, UInt_t centC=k10) const
Bool_t isNewFile(TString &name)
Definition: hloop.cc:1013
Int_t getSelectedParticleCandMult(Int_t s) const
Float_t getEventPlaneParameter(UInt_t flag=kQx, Bool_t corr=kFALSE)
Bool_t isFlagSet(UInt_t flag, UInt_t status)
Int_t getCentralityClassNbins(UInt_t centC=k10)
Bool_t goodSector(Int_t sector)
Definition: hloop.h:188
TH2D * makeEPresolution(TProfile2D *hist, Bool_t calcChi=kFALSE)
vector< Float_t > vQPhi
Float_t getFWSumZ(Int_t minZ=1, Int_t maxZ=99, UInt_t SubEvent=0)
Short_t getSystemUsed() const
void getXYZLab(Float_t &x, Float_t &y, Float_t &z)
Definition: hwallhit.h:42
Float_t getTime(void) const
Definition: hwallhit.h:35
Float_t getSmearValue(Int_t size=1)
Float_t phi
Definition: drawAccCuts.C:15
Int_t getSumPrimaryParticleCandMult() const
TString getStringFWCutValues(UInt_t e)
const Cat_t catParticleEvtInfo
Definition: hparticledef.h:10
Bool_t loadEventPlaneCorrectionHist()
vector< TH1 * > fCentralityPercentileHist
vector< vector< TH1 * > > fFWCutValuesHist
Int_t getPID() const
Bool_t fillQVectors(UInt_t EPcorr=kNoCorrection, UInt_t nHarmonic=1)
Bool_t addEventPlaneCorrectionHist(TProfile2D *hist, UInt_t epParam=kQx)
UInt_t getTBit(void)
Definition: heventheader.h:136
vector< vector< TH1 * > > fCentralityHist
Float_t getCharge(void) const
Definition: hwallhit.h:36
Bool_t printCentralityClass(TString estimator)
Float_t getEventPlaneWeight(UInt_t EPcorr=kNoCorrection, UInt_t SubEvent=0, UInt_t nHarmonic=1)
#define MAXFWCELLS
Double_t ModifiedBesselI(Int_t n, Double_t x) const
TH1F * getEventPlaneCorrectionHist(UInt_t flag) const
TString getStringEventPlaneParameter(UInt_t e)
Int_t getCentralityClass10(TString estimator)
vector< Int_t > iFWHitvector
Int_t getSumTofMult() const
Int_t getSumTofMultCut() const
Int_t getSumRpcMultHitCut() const
const Cat_t catParticleCand
Definition: hparticledef.h:9
vector< SimpleQVector * > arrayOfHits
Bool_t isFlagAND(Int_t num,...)