Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
Femto1DCFAnaMapPairsDumped.cxx
1/*
2 * Femto1DCFAnaMapPairs.cxx
3 *
4 * Created on: 27 lip 2019
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "Femto1DCFAnaMapPairsDumped.h"
11
12#include "DividedHisto.h"
13#include "Femto1DCFAnaMapMC.h"
14#include "Femto1DMapGenerator.h"
15#include "FemtoFreezoutGenerator.h"
16#include "FemtoMiniPair.h"
17#include "FemtoPair.h"
18#include "FemtoSourceModel.h"
19#include "FemtoWeightGenerator.h"
20
21#include <TBranch.h>
22#include <TClonesArray.h>
23#include <TCollection.h>
24#include <TDirectoryFile.h>
25#include <TFile.h>
26#include <TH2.h>
27#include <TKey.h>
28#include <TList.h>
29#include <TMathBase.h>
30#include <TNamed.h>
31#include <TObjArray.h>
32#include <TObject.h>
33#include <TTree.h>
34
35namespace Hal {
36 Femto1DCFAnaMapPairsDumped::Femto1DCFAnaMapPairsDumped() :
37 Femto1DCFAnaMapMC(kFALSE),
38 fInFile("pairs.root"),
39 fFile(nullptr),
40 fTree(nullptr),
41 fPairsSignal(nullptr),
42 fPairsBackground(nullptr),
43 fUseBackground(kFALSE),
44 fWeightedBackround(kFALSE),
45 fMaxEvents(-1) {}
46
47 void Femto1DCFAnaMapPairsDumped::SetInFile(TString filename) { fInFile = filename; }
48
49 void Femto1DCFAnaMapPairsDumped::Run(Int_t start, Int_t end) {
50 if (start == end) {
51 start = 0;
52 end = fTree->GetEntries();
53 }
54 for (int i = start; i < end; i++) {
55 fTree->GetEntry(i);
56 if (!fUseBackground) {
57 MakePairsWeighted(fPairsSignal, kTRUE);
58 MakePairsUnWeighted(fPairsSignal, kFALSE);
59 } else {
60 MakePairsWeighted(fPairsSignal, kTRUE);
61 if (fWeightedBackround) {
62 MakePairsWeighted(fPairsBackground, kFALSE);
63 } else {
64 MakePairsUnWeighted(fPairsBackground, kFALSE);
65 }
66 }
67 }
68 }
69
70 Bool_t Femto1DCFAnaMapPairsDumped::Init() {
71 fFile = new TFile(fInFile);
72 if (fFile->IsZombie()) return kFALSE;
73 TList* list = fFile->GetListOfKeys();
74 for (int i = 0; i < list->GetEntries(); i++) {
75 TKey* key = (TKey*) list->At(i);
76 TObject* obj = fFile->Get(key->GetName());
77 if (obj->InheritsFrom("TTree")) {
78 fTree = (TTree*) obj;
79 break;
80 }
81 }
82 if (fTree == nullptr) return kFALSE;
83 TBranch* br = fTree->GetBranch("FemtoSignal.");
84 fPairsSignal = new TClonesArray("Hal::FemtoMiniPair");
85 br->SetAddress(&fPairsSignal);
86 if (fUseBackground) {
87 br = fTree->GetBranch("FemtoBackground.");
88 fPairsBackground = new TClonesArray("Hal::FemtoMiniPair");
89 br->SetAddress(&fPairsBackground);
90 }
91 return Femto1DCFAnaMapMC::Init();
92 }
93
94 void Femto1DCFAnaMapPairsDumped::MakePairsWeighted(TClonesArray* arr, Bool_t num) {
95 TH2* hist;
96 if (num) {
97 hist = (TH2*) fMap->GetNum();
98 } else {
99 hist = (TH2*) fMap->GetDen();
100 }
101
102 for (int i = 0; i < arr->GetEntriesFast(); i++) {
103 FemtoMiniPair* mini = (FemtoMiniPair*) arr->UncheckedAt(i);
104 *fPair = *mini;
105 fPair->Compute();
106 Double_t k = fPair->GetT();
107 if (fIgnoreSign) { k = TMath::Abs(k); }
108 for (int r_bin = 0; r_bin < fRBins; r_bin++) {
109 Double_t R = fRadiiBins[r_bin];
110 fGenerator->GetSourceModel()->SetRadius(R);
111 fGenerator->GenerateFreezoutCooordinates(fPair);
112 Double_t weight = fWeight->GenerateWeight(fPair);
113 hist->Fill(k, R, weight);
114 }
115 }
116 }
117
118 void Femto1DCFAnaMapPairsDumped::MakePairsUnWeighted(TClonesArray* arr, Bool_t num) {
119 TH2* hist;
120 if (num) {
121 hist = (TH2*) fMap->GetNum();
122 } else {
123 hist = (TH2*) fMap->GetDen();
124 }
125 for (int i = 0; i < arr->GetEntriesFast(); i++) {
126 FemtoMiniPair* mini = (FemtoMiniPair*) arr->UncheckedAt(i);
127 *fPair = *mini;
128 fPair->Compute();
129 Double_t k = fPair->GetT();
130 if (fIgnoreSign) { k = TMath::Abs(k); }
131 for (double r_bin = 0; r_bin < fRBins; r_bin++) {
132 Double_t R = fRadiiBins[r_bin];
133 hist->Fill(k, R, 1);
134 }
135 }
136 }
137
138 Femto1DCFAnaMapPairsDumped::~Femto1DCFAnaMapPairsDumped() {
139 // TODO Auto-generated destructor stub
140 }
141} // namespace Hal