Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
Femto3DCFQinv.cxx
1/*
2 * Femto3DCFQinv.cxx
3 *
4 * Created on: 19 lut 2019
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "Femto3DCFQinv.h"
11#include "FemtoConst.h"
12#include "FemtoPair.h"
13#include <TBrowser.h>
14#include <TCollection.h>
15#include <TH3.h>
16#include <TMathBase.h>
17#include <TString.h>
18#include <TVirtualPad.h>
19
20namespace Hal {
21 Femto3DCFQinv::Femto3DCFQinv(TString name, Femto::EKinematics frame) : Femto3DCF(name, frame), fPureNum(NULL), fQinvNum(NULL) {
22 if (GetNum()) {
23 fPureNum = (TH3D*) GetNum()->Clone("qinvN");
24 fQinvNum = (TH3D*) GetNum()->Clone("qinvD");
25 }
26 }
27
28 Femto3DCFQinv::Femto3DCFQinv(TString name,
29 Int_t binsX,
30 Double_t minX,
31 Double_t maxX,
32 Int_t binsY,
33 Double_t minY,
34 Double_t maxY,
35 Int_t binsZ,
36 Double_t minZ,
37 Double_t maxZ,
38 Femto::EKinematics frame) :
39 Femto3DCF(name, binsX, minX, maxX, binsY, minY, maxY, binsZ, minZ, maxZ, frame) {
40 fPureNum = (TH3D*) GetNum()->Clone("qinvN");
41 fQinvNum = (TH3D*) GetNum()->Clone("qinvD");
42 }
43
44 Femto3DCFQinv::Femto3DCFQinv(const Femto3DCFQinv& other) : Femto3DCF(other) {
45 fPureNum = (TH3D*) other.fPureNum->Clone("qinvN");
46 fQinvNum = (TH3D*) other.fQinvNum->Clone("qinvD");
47 }
48
49 void Femto3DCFQinv::FillNumObj(TObject* obj) {
50 FemtoPair* pair = (FemtoPair*) obj;
51 Double_t x = pair->GetX();
52 Double_t y = pair->GetY();
53 Double_t z = pair->GetZ();
54 if (pair->IsAbs()) {
55 x = TMath::Abs(x);
56 y = TMath::Abs(y);
57 z = TMath::Abs(z);
58 }
59 ((TH3*) fNum)->Fill(x, y, z, pair->GetWeight());
60 ((TH3*) fPureNum)->Fill(x, y, z);
61 ((TH3*) fQinvNum)->Fill(x, y, z, pair->GetT());
62 }
63
64 TH3D* Femto3DCFQinv::GetQinvHist(TString opt) const {
65 if (opt == "av") {
66 TH3D* div = (TH3D*) fQinvNum->Clone("ratio");
67 div->Divide(fPureNum);
68 return div;
69 }
70 return NULL;
71 }
72
73 Long64_t Femto3DCFQinv::Merge(TCollection* collection) {
74 if (collection) {
75 Femto3DCFQinv* pack = NULL;
76 TIter iterator(collection);
77 while ((pack = (Femto3DCFQinv*) iterator())) {
78 Add(pack);
79 }
80 }
81 return 1;
82 }
83
84 Femto3DCFQinv& Femto3DCFQinv::operator=(const Femto3DCFQinv& other) {
85 if (&other != this) {
86 DividedHisto3D::operator=(other);
87 if (fPureNum) delete fPureNum;
88 if (fQinvNum) delete fQinvNum;
89 fPureNum = (TH3D*) other.fPureNum->Clone();
90 fQinvNum = (TH3D*) other.fQinvNum->Clone();
91 }
92 return *this;
93 }
94
95 void Femto3DCFQinv::Add(const Object* h) {
96 Femto3DCF::Add(h);
97 const Femto3DCFQinv* cf = dynamic_cast<const Femto3DCFQinv*>(h);
98 if (cf) {
99 fPureNum->Add(cf->fPureNum);
100 fQinvNum->Add(cf->fQinvNum);
101 }
102 }
103
104 void Femto3DCFQinv::AddNum(TH1* num, Option_t* opt) {
105 Femto3DCF::AddNum(num, opt);
106 if (fPureNum) {
107 delete fPureNum;
108 delete fQinvNum;
109 }
110 fPureNum = (TH3D*) GetNum()->Clone("qinvN");
111 fQinvNum = (TH3D*) GetNum()->Clone("qinvD");
112 }
113
114 void Femto3DCFQinv::AddDen(TH1* den, Option_t* opt) { Femto3DCF::AddDen(den, opt); }
115
116 void Femto3DCFQinv::Browse(TBrowser* b) {
117 TVirtualPad* c1 = gPad;
118 gPad->Clear();
119 Draw("all");
120 gPad = c1;
121 b->SetDrawOption("colz");
122 b->Add(fNum);
123 b->Add(fDen);
124 b->Add(GetQinvHist("av"));
125 }
126
127 Femto3DCFQinv::Femto3DCFQinv() : Femto3DCF(), fPureNum(nullptr), fQinvNum(nullptr) {}
128
129 Femto3DCFQinv::~Femto3DCFQinv() {
130 if (fQinvNum) {
131 delete fQinvNum;
132 delete fPureNum;
133 }
134 }
135
136 void Femto3DCFQinv::AddScaled(const Hal::DividedHisto1D& other, Double_t scale) {
138 auto other2 = dynamic_cast<const Hal::Femto3DCFQinv&>(other);
139 fPureNum->Add(other2.fPureNum, scale);
140 fQinvNum->Add(other2.fQinvNum, scale);
141 }
142
143} // namespace Hal
virtual void AddScaled(const DividedHisto1D &other, Double_t scale=1)
virtual void AddScaled(const Hal::DividedHisto1D &other, Double_t scale)
Double_t GetT() const
Definition FemtoPair.h:322
Double_t GetY() const
Definition FemtoPair.h:312
Double_t GetWeight() const
Definition FemtoPair.h:282
Double_t GetZ() const
Definition FemtoPair.h:317
Double_t GetX() const
Definition FemtoPair.h:307