Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoDumpPairAna.cxx
1/*
2 * FemtoDumpPairAna.cxx
3 *
4 * Created on: 9 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 "FemtoDumpPairAna.h"
11
12
13#include "Cout.h"
14#include "DataManager.h"
15#include "DividedHisto.h"
16#include "FemtoCorrFunc.h"
17#include "FemtoPair.h"
18#include "MemoryMapManager.h"
19#include "ObjectMatrix.h"
20#include "TwoTrackAna.h"
21
22#include <TAxis.h>
23#include <TClonesArray.h>
24#include <TFile.h>
25#include <TH1.h>
26#include <TObjArray.h>
27
28
29namespace Hal {
30 FemtoDumpPairAna::FemtoDumpPairAna() : fBinLimit(10000), fLimitsN(), fLimitsD(), fWriteBackground(kFALSE) {}
31
32 void FemtoDumpPairAna::ProcessFemtoPair() {
33 fFemtoPair->Compute();
34 Int_t bin = fGrouping.GetBin(fFemtoPair);
35 if (bin < 0) return;
36 if (fLimitsN[bin] < fBinLimit) {
37 FemtoMicroPair* mini = (FemtoMicroPair*) fSignalPairs[bin]->ConstructedAt(fSignalPairs[bin]->GetEntriesFast());
38 *mini = *fFemtoPair;
39 fLimitsN[bin]++;
40 }
41 }
42
43 void FemtoDumpPairAna::ProcessFemtoPair_Rotated() {
44 fFemtoPair->Compute_Rotated();
45 Int_t bin = fGrouping.GetBin(fFemtoPair);
46 if (bin < 0) return;
47 if (fLimitsN[bin] < fBinLimit) {
48 FemtoMicroPair* mini = (FemtoMicroPair*) fBackgroundPairs[bin]->ConstructedAt(fBackgroundPairs[bin]->GetEntriesFast());
49 *mini = *fFemtoPair;
50 fLimitsD[bin]++;
51 }
52 }
53
54 void FemtoDumpPairAna::ProcessFemtoPair_Hemisphere() {
55 fFemtoPair->Compute_Hemisphere();
56 Int_t bin = fGrouping.GetBin(fFemtoPair);
57 if (bin < 0) return;
58 if (fLimitsN[bin] < fBinLimit) {
59 FemtoMicroPair* mini = (FemtoMicroPair*) fBackgroundPairs[bin]->ConstructedAt(fBackgroundPairs[bin]->GetEntriesFast());
60 *mini = *fFemtoPair;
61 fLimitsD[bin]++;
62 }
63 }
64
65 void FemtoDumpPairAna::ProcessFemtoPair_Mixed() {
66 fFemtoPair->Compute();
67 Int_t bin = fGrouping.GetBin(fFemtoPair);
68 if (bin < 0) return;
69 if (fLimitsN[bin] < fBinLimit) {
70 FemtoMicroPair* mini = (FemtoMicroPair*) fBackgroundPairs[bin]->ConstructedAt(fBackgroundPairs[bin]->GetEntriesFast());
71 *mini = *fFemtoPair;
72 fLimitsD[bin]++;
73 }
74 }
75
76 Task::EInitFlag FemtoDumpPairAna::Init() {
77 Task::EInitFlag stat = FemtoBasicAna::Init();
78 if (stat == Task::EInitFlag::kFATAL) return stat;
79 DataManager* mngr = DataManager::Instance();
80 FemtoCorrFunc* f = (FemtoCorrFunc*) fCFs->At(0, 0);
81 DividedHisto1D* dummy = f->GetCF(0);
82 Int_t bins;
83 Double_t min, max;
84
85 // grouping enabled
86 if (dummy->GetNum()->InheritsFrom("TH3")) { // 3d histo group by q-long/ k-long
87 fGrouping.SetGroupByKLong();
88 Std::GetAxisPar(*dummy->GetNum(), bins, min, max, "z");
89 } else { // 1d histo group by k* qinv
90 fGrouping.GroupByKStar();
91 Std::GetAxisPar(*dummy->GetNum(), bins, min, max, "x");
92 }
93 fGrouping.SetFrame(fFemtoPair->GetFrame());
94 fGrouping.SetAxis(bins, min, max);
95 fLimitsN.MakeBigger(bins);
96 fLimitsD.MakeBigger(bins);
97 auto vec = fGrouping.GetBranchesByValue(0, 0, kTRUE); // 0,0 -> get all branches
98 int idx = 0;
99 for (auto branchName : vec) {
100 fSignalPairs.push_back(new TClonesArray("Hal::FemtoMicroPair", 100));
101 mngr->Register(branchName, "FemtoPairs", fSignalPairs[idx++], kTRUE);
102 }
103 if (fWriteBackground) {
104 idx = 0;
105 vec = fGrouping.GetBranchesByValue(0, 0, kFALSE);
106 for (auto branchName : vec) {
107 fBackgroundPairs.push_back(new TClonesArray("Hal::FemtoMicroPair", 100));
108 mngr->Register(branchName, "FemtoPairs", fBackgroundPairs[idx++], kTRUE);
109 }
110 }
111
112 if (IdenticalParticles())
113 fBackgroundMode = kNoBackgroundID;
114 else
115 fBackgroundMode = kNoBackgroundNID;
116 return Task::EInitFlag::kSUCCESS;
117 }
118
119 FemtoDumpPairAna::~FemtoDumpPairAna() {}
120
121 void FemtoDumpPairAna::Exec(Option_t* opt) {
122 for (auto i : fSignalPairs)
123 i->Clear();
124
125 for (auto i : fBackgroundPairs)
126 i->Clear();
127
128 FemtoBasicAna::Exec(opt);
129 }
130
131 void FemtoDumpPairAna::FinishTask() {
133 // Package* pack = Report();
134 GoToDir("HalInfo");
135 fGrouping.Clone()->Write();
136 gFile->cd();
137 }
138
139 //==================================================================== group class
140} // namespace Hal
void Register(const char *name, const char *folderName, TNamed *obj, Bool_t toFile)
virtual void FinishTask()
DividedHisto1D * GetCF(Int_t i) const