Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
Femto1DMapGenerator.cxx
1/*
2 * Femto1DMapGenerator.cxx
3 *
4 * Created on: 26 sie 2021
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#include "Femto1DMapGenerator.h"
10
11#include "CorrFitMapKstarRstar.h"
12#include "Cout.h"
13#include "Femto1DCF.h"
14#include "FemtoFreezoutGenerator.h"
15#include "FemtoPair.h"
16#include "FemtoWeightGenerator.h"
17#include "Std.h"
18
19#include <TDatabasePDG.h>
20#include <TFile.h>
21#include <iostream>
22
23namespace Hal {
24 Femto1DMapGenerator::Femto1DMapGenerator(Bool_t fake) :
25 fRMin(0.5),
26 fRMax(10.5),
27 fKStarMin(0),
28 fKStarMax(1.0),
29 fMass12(-1),
30 fMass22(-1),
31 fRBins(10),
32 fKStarBins(100),
33 fPid1(0),
34 fPid2(0),
35 fIgnoreSign(kFALSE),
36 fUseFake(fake),
37 fKinematics(Femto::EKinematics::kPRF),
38 fPair(nullptr),
39 fGenerator(nullptr),
40 fWeight(nullptr),
41 fMap(nullptr) {}
42
43 Femto1DMapGenerator::Femto1DMapGenerator() : Femto1DMapGenerator(kFALSE) {}
44
45 void Femto1DMapGenerator::SetCorrFctn(const Femto1DCF& cf) {
46 fKinematics = cf.GetFrame();
47 fKStarBins = cf.GetNum()->GetXaxis()->GetNbins();
48 fKStarMin = cf.GetNum()->GetXaxis()->GetBinLowEdge(1);
49 fKStarMax = cf.GetNum()->GetXaxis()->GetBinUpEdge(fKStarBins);
50 }
51
52 void Femto1DMapGenerator::SetPid(Int_t pid1, Int_t pid2) {
53 TDatabasePDG* pid = TDatabasePDG::Instance();
54 TParticlePDG* p1 = pid->GetParticle(pid1);
55 TParticlePDG* p2 = pid->GetParticle(pid2);
56 if (p1 == nullptr) {
57 Cout::PrintInfo(Form("Femto1DMapGenerator::SetPid cannot find PID1 = %i", pid1), EInfo::kLowWarning);
58 return;
59 };
60 if (p2 == nullptr) {
61 Cout::PrintInfo(Form("Femto1DMapGenerator::SetPid cannot find PID1 = %i", pid2), EInfo::kLowWarning);
62 return;
63 };
64 fPid1 = pid1;
65 fPid2 = pid2;
66 fMass12 = p1->Mass() * p1->Mass();
67 fMass22 = p2->Mass() * p2->Mass();
68 }
69
70 void Femto1DMapGenerator::SetRBins(Int_t bins, Double_t rmin, Double_t rmax, Bool_t center) {
71 if (center) {
72 fRBins = bins + 1;
73 fRadiiBins.MakeBigger(bins + 1);
74 Double_t dr = (rmax - rmin) / (Double_t)(bins);
75 fRMin = rmin - dr * 0.5;
76 fRMax = rmax + dr * 0.5;
77 for (int i = 0; i < fRBins; i++) {
78 fRadiiBins[i] = fRMin + dr * i + dr * 0.5;
79 }
80 } else {
81 fRBins = bins;
82 fRadiiBins.MakeBigger(bins);
83 fRMin = rmin;
84 fRMax = rmax;
85 Double_t dr = (fRMax - fRMin) / (Double_t)(bins);
86 for (int i = 0; i < fRBins; i++) {
87 fRadiiBins[i] = fRMin + dr * i + dr * 0.5;
88 }
89 }
90 }
91
92 void Femto1DMapGenerator::SetGenerator(FemtoFreezoutGenerator& gen) {
93 if (fGenerator) delete fGenerator;
94 fGenerator = gen.MakeCopy();
95 }
96
97 void Femto1DMapGenerator::SetWeightGenerator(FemtoWeightGenerator& calc) {
98 if (fWeight) delete fWeight;
99 fWeight = calc.MakeCopy();
100 }
101
102 Bool_t Femto1DMapGenerator::Init() {
103 fPair = Femto::MakePair(fKinematics, fUseFake);
104 fMap = new DividedHisto2D("map", fKStarBins, fKStarMin, fKStarMax, fRBins, fRMin, fRMax, 'D');
105 if (fWeight == nullptr) {
106 Cout::PrintInfo("Femto1DMapGenerator::Init lack of weight", EInfo::kCriticalError);
107 return kFALSE;
108 }
109 if (fGenerator == nullptr) {
110 Cout::PrintInfo("Femto1DMapGenerator::Init lack of generator", EInfo::kCriticalError);
111 return kFALSE;
112 }
113 if (fRadiiBins.GetSize() < 1) {
114 Cout::PrintInfo("Femto1DMapGenerator::Init lack of rbins ! did you call SetRBins?", EInfo::kCriticalError);
115 return kFALSE;
116 }
117 fPair->SetPdg1(fPid1);
118 fPair->SetPdg2(fPid2);
119 if (fIgnoreSign) fPair->UseAbs();
120 return kTRUE;
121 }
122
123 void Femto1DMapGenerator::SaveMap(TString filename) {
124 TFile* f = new TFile(filename, "recreate");
125 CorrFitMapKstarRstarDiv* map = new CorrFitMapKstarRstarDiv(*fMap, fKinematics);
126 map->Write("map");
127 f->Close();
128 }
129
130 Femto1DMapGenerator::~Femto1DMapGenerator() {
131 if (fPair) delete fPair;
132 if (fGenerator) delete fGenerator;
133 if (fWeight) delete fWeight;
134 if (fMap) delete fMap;
135 }
136} // namespace Hal
Femto::EKinematics GetFrame() const
Definition Femto1DCF.h:66
virtual FemtoFreezoutGenerator * MakeCopy() const =0
virtual FemtoWeightGenerator * MakeCopy() const