Project
Loading...
Searching...
No Matches
FlowMapper.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13#include "TH1D.h"
14#include "TH3D.h"
15#include "TF1.h"
16#include <fairlogger/Logger.h>
17
18namespace o2
19{
20namespace eventgen
21{
22
23/*****************************************************************/
24/*****************************************************************/
25
27{
28 // base constructor. Creates cumulative function only so that's already in place but not much else.
29 mCumulative = std::make_unique<TF1>("mCumulative", "x+[0]*TMath::Sin(2*x)", 0, 2 * TMath::Pi());
30 mBinsPhi = 4000; // a first guess
31}
32
33void FlowMapper::CreateLUT(TH1D* mhv2vsPt, TH1D* mhEccVsB)
34{
35 if (!mhv2vsPt) {
36 LOG(fatal) << "You did not specify a valid v2 vs pT histogram!";
37 return;
38 }
39 if (!mhEccVsB) {
40 LOG(fatal) << "You did not specify a valid ecc vs B histogram!";
41 return;
42 }
43 LOG(info) << "Proceeding to creating a look-up table...";
44 const Long_t nbinsB = mhEccVsB->GetNbinsX();
45 const Long_t nbinsPt = mhv2vsPt->GetNbinsX();
46 const Long_t nbinsPhi = mBinsPhi; // constant in this context necessary
47 std::vector<double> binsB(nbinsB + 1, 0);
48 std::vector<double> binsPt(nbinsPt + 1, 0);
49 std::vector<double> binsPhi(nbinsPhi + 1, 0);
50
51 for (int ii = 0; ii < nbinsB + 1; ii++) {
52 binsB[ii] = mhEccVsB->GetBinLowEdge(ii + 1);
53 }
54 for (int ii = 0; ii < nbinsPt + 1; ii++) {
55 binsPt[ii] = mhv2vsPt->GetBinLowEdge(ii + 1);
56 }
57 for (int ii = 0; ii < nbinsPhi + 1; ii++) {
58 binsPhi[ii] = static_cast<Double_t>(ii) * 2 * TMath::Pi() / static_cast<Double_t>(nbinsPhi);
59 }
60
61 // std::make_unique<TH1F>("hSign", "Sign of electric charge;charge sign", 3, -1.5, 1.5);
62
63 mhLUT = std::make_unique<TH3D>("mhLUT", "", nbinsB, binsB.data(), nbinsPt, binsPt.data(), nbinsPhi, binsPhi.data());
64 mhLUT->SetDirectory(nullptr); // just in case context is incorrect
65
66 // loop over each centrality (b) bin
67 for (int ic = 0; ic < nbinsB; ic++) {
68 // loop over each pt bin
69 for (int ip = 0; ip < nbinsPt; ip++) {
70 // find target v2 value and set cumulative for inversion
71 double v2target = mhv2vsPt->GetBinContent(ip + 1) * mhEccVsB->GetBinContent(ic + 1);
72 LOG(info) << "At b ~ " << mhEccVsB->GetBinCenter(ic + 1) << ", pt ~ " << mhv2vsPt->GetBinCenter(ip + 1) << ", ref v2 is " << mhv2vsPt->GetBinContent(ip + 1) << ", scale is " << mhEccVsB->GetBinContent(ic + 1) << ", target v2 is " << v2target << ", inverting...";
73 mCumulative->SetParameter(0, v2target); // set up
74 for (Int_t ia = 0; ia < nbinsPhi; ia++) {
75 // Look systematically for the X value that gives this Y
76 // There are probably better ways of doing this, but OK
77 Double_t lY = mhLUT->GetZaxis()->GetBinCenter(ia + 1);
78 Double_t lX = lY; // a first reasonable guess
79 Bool_t lConverged = kFALSE;
80 while (!lConverged) {
81 Double_t lDistance = mCumulative->Eval(lX) - lY;
82 if (TMath::Abs(lDistance) < mPrecision) {
83 lConverged = kTRUE;
84 break;
85 }
86 Double_t lDerivativeValue = mDerivative / (mCumulative->Eval(lX + mDerivative) - mCumulative->Eval(lX));
87 lX = lX - lDistance * lDerivativeValue * 0.25; // 0.5: speed factor, don't overshoot but control reasonable
88 }
89 mhLUT->SetBinContent(ic + 1, ip + 1, ia + 1, lX);
90 }
91 }
92 }
93}
94
95Double_t FlowMapper::MapPhi(Double_t lPhiInput, Double_t b, Double_t pt)
96{
97 Int_t lLowestPeriod = TMath::Floor(lPhiInput / (2 * TMath::Pi()));
98 Double_t lPhiOld = lPhiInput - 2 * lLowestPeriod * TMath::Pi();
99 Double_t lPhiNew = lPhiOld;
100
101 // Avoid interpolation problems in dimension: pT
102 Double_t lMaxPt = mhLUT->GetYaxis()->GetBinCenter(mhLUT->GetYaxis()->GetNbins());
103 Double_t lMinPt = mhLUT->GetYaxis()->GetBinCenter(1);
104 if (pt > lMaxPt) {
105 pt = lMaxPt; // avoid interpolation problems at edge
106 }
107
108 Double_t phiWidth = mhLUT->GetZaxis()->GetBinWidth(1); // any bin, assume constant
109
110 // Valid if not at edges. If at edges, that's ok, do not map
111 bool validPhi = lPhiNew > phiWidth / 2.0f && lPhiNew < 2.0 * TMath::Pi() - phiWidth / 2.0f;
112
113 // If at very high b, do not map
114 bool validB = b < mhLUT->GetXaxis()->GetBinCenter(mhLUT->GetXaxis()->GetNbins());
115 Double_t minB = mhLUT->GetXaxis()->GetBinCenter(1);
116
117 if (validPhi && validB) {
118
119 Double_t scaleFactor = 1.0; // no need if not special conditions
120 if (pt < lMinPt) {
121 scaleFactor *= pt / lMinPt; // downscale the difference, zero at zero pT
122 pt = lMinPt;
123 }
124 if (b < minB) {
125 scaleFactor *= b / minB; // downscale the difference, zero at zero b
126 b = minB;
127 }
128 lPhiNew = mhLUT->Interpolate(b, pt, lPhiOld);
129
130 lPhiNew = scaleFactor * lPhiNew + (1.0 - scaleFactor) * lPhiOld;
131 }
132 return lPhiNew + 2.0 * lLowestPeriod * TMath::Pi();
133}
134
135/*****************************************************************/
136/*****************************************************************/
137
138} /* namespace eventgen */
139} /* namespace o2 */
140
ClassImp(o2::eventgen::FlowMapper)
std::unique_ptr< TF1 > mCumulative
Definition FlowMapper.h:61
std::unique_ptr< TH3D > mhLUT
Definition FlowMapper.h:64
void CreateLUT(TH1D *mhv2vsPt, TH1D *mhEccVsB)
Double_t MapPhi(Double_t lPhiInput, Double_t b, Double_t pt)
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"