Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoSHCF.h
1/*
2 * FemtoSHCF.h
3 *
4 * Created on: 28 kwi 2016
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#ifndef HALFEMTOSHCF_H_
10#define HALFEMTOSHCF_H_
11
12#include "FemtoConst.h"
13
14#include <initializer_list>
15
16#ifdef __CIA__
17#include "CorrFitMaskSH.h"
18#endif
19#include "Array.h"
20#include "DividedHisto.h"
21#include "FemtoYlmIndexes.h"
22#include "FemtoYlmMath.h"
23#include "Style.h"
24
25#include <TMath.h>
26#include <TString.h>
27#include <complex>
28
33namespace Hal {
34#ifdef __CIA__
35 class CorrFitSHCF;
36#endif
37 class FemtoSHSlice;
38 class FemtoYlmSolver;
39 class FemtoCFPainter;
40
44 class FemtoSHCF : public DividedHisto1D {
45#ifdef __CIA__
46 friend class CorrFitSCHF;
47#endif
49 friend class FemtoSHSlice;
50 friend class FemtoYlmSolver;
51 friend class FemtoSHCFPainter;
52 const Int_t fMaxJM;
53 TH1D** fNumReal = {nullptr}; //[fMaxJM] Real parts of Ylm components of the numerator
54 TH1D** fNumImag = {nullptr}; // [fMaxJM]Imaginary parts of Ylm components of the numerator
55 TH1D** fDenReal = {nullptr}; //[fMaxJM] Real parts of Ylm components of the denominator
56 TH1D** fDenImag = {nullptr}; //[fMaxJM] Imaginary parts of Ylm components of the denominator
57 Femto::EKinematics fFrame;
58
59 TH1D** fCFReal = {nullptr}; //[fMaxJM] real CF's
60 TH1D** fCFImag = {nullptr}; //[fMaxJM] img CF's
61 Array_3<Double_t> fCovNum;
62 Array_3<Double_t> fCovDen;
63 Array_3<Double_t> fCovCf;
64 Double_t fNormPurity; //
65 Double_t fNormRadius; //
66 Double_t fNormBohr; //
67 TH3D* fCfcov = {nullptr}; //
68 FemtoCFPainter* fPainter = {nullptr};
69 FemtoYlmIndexes fLmVals;
70 FemtoYlmMath fLmMath;
71 Bool_t fColzSet = {kFALSE};
72 Color_t fColRe = {kBlue};
73 Color_t fColIm = {kRed};
74
75 Double_t Sqr(Double_t val1, Double_t val2) const;
76 TH1D* Histo(int ilm, int em, Option_t* opt, Double_t scale) const;
77 void PackCfcCovariance();
78 void UnpackCovariances() {};
87 TH1D* GetHisto(int el, int em, Bool_t norm, Option_t* opt = "") const;
88
89 protected:
94 Int_t GetMaxJM() const { return fMaxJM; };
99 virtual void FastAdd(const FemtoSHCF* obj);
100 virtual void ApplyStyle(const Hal::HistoStyle& h);
101
102 public:
105 FemtoSHCF();
115 FemtoSHCF(TString name, Int_t maxL, Int_t bins, Double_t min, Double_t max, Femto::EKinematics kinematics);
120 FemtoSHCF(const FemtoSHCF& other);
121 void SetCFColz(Color_t re, Color_t im) {
122 fColRe = re;
123 fColIm = im;
124 fColzSet = kTRUE;
125 }
126 virtual void FillNumObj(TObject* obj);
127 virtual void FillDenObj(TObject* obj);
128 FemtoYlmIndexes GetYlmIndexes() const { return fLmVals; }
134 /* TODO implenet this function in proper way
135 FemtoSHCF(Femto3DCF *cf, Int_t maxL=3);
136 */
141 void SetComment(TString comment) { fComment = comment; };
146 void SetNormPurity(Double_t purity) { fNormPurity = purity; }
151 void SetNormRadius(Double_t radius) { fNormRadius = radius; };
156 void SettNormBohr(Double_t bohr) { fNormBohr = bohr; };
164 void AddLabel(TString label);
172 void RecalculateCF(Int_t debugBin = -1, Bool_t suwm = kFALSE);
179 void SetNumRe(TH1D** histograms, Bool_t clone = kTRUE);
186 void SetNumIm(TH1D** histograms, Bool_t clone = kTRUE);
193 void SetDenRe(TH1D** histograms, Bool_t clone = kTRUE);
200 void SetDenIm(TH1D** histograms, Bool_t clone = kTRUE);
205 void Draw(Option_t* opt = "");
213 void NormalizeBy(int el, int em, Option_t* opt = "re");
218 void SetScale(Double_t scale) { fScale = scale; };
223 void Browse(TBrowser* b);
224#ifdef __CIA__
228 void Fit(CorrFitSHCF* fit);
234 void FitDummy(CorrFitSHCF* fit);
235#endif
241 TH1D* GetRawCFRe(int pos) const { return fCFReal[pos]; };
247 TH1D* GetRawCFIm(int pos) const { return fCFImag[pos]; };
254 TH1D* GetCFRe(int el, int em) const;
261 TH1D* GetCFIm(int el, int em) const;
268 TH1D* GetNumRe(int el, int em) const;
275 TH1D* GetNumIm(int el, int em) const;
282 TH1D* GetDenRe(int el, int em) const;
289 TH1D* GetDenIm(int el, int em) const;
294 TH3D* GetCovCF() const;
300 TH1D* GetCubic(Option_t* opt, Bool_t R = kFALSE) const;
301 Femto::EKinematics GetFrame() const { return fFrame; }
306 Int_t GetLMax() const { return TMath::Sqrt(fMaxJM - 1); };
307 void Rebin(Int_t ngroup, Option_t* opt = "");
308 TH2D* GetMaskDraw() const;
309#ifdef __CIA__
314 CorrFitMaskSH MakeEmptyMask() const;
315#endif
316 Array_1<Float_t>* ExportToFlatNum() const;
317 void ExportIntoToFlatNum(Array_1<Float_t>* output) const;
318 void ExportIntoToFlatNumValkyria(Array_1<Float_t>* output) const;
319#ifdef __CIA__
326#endif
327 virtual void ImportSlice(Array_1<Float_t>* array, Int_t bin);
328 virtual void Add(const Object* pack);
329 virtual Long64_t Merge(TCollection* collection);
330 void MakeDummyCov();
331 Array_3<Double_t> GetCovNum() { return fCovNum; }
332 virtual TString HTMLExtract(Int_t counter = 0, TString dir = " ") const;
333 virtual TObject* GetSpecial(TString opt) const;
338 FemtoCFPainter* GetPainter() const { return fPainter; } // KURWA
339 virtual ~FemtoSHCF();
340 ClassDef(FemtoSHCF, 6)
341 };
342} // namespace Hal
343#endif /* HALFEMTOSHCF_H_ */
TH1D * GetCFRe(int el, int em) const
void NormalizeBy(int el, int em, Option_t *opt="re")
void SetDenIm(TH1D **histograms, Bool_t clone=kTRUE)
virtual void FillNumObj(TObject *obj)
virtual void ApplyStyle(const Hal::HistoStyle &h)
TH1D * GetCubic(Option_t *opt, Bool_t R=kFALSE) const
FemtoCFPainter * GetPainter() const
Definition FemtoSHCF.h:338
TH1D * GetRawCFRe(int pos) const
Definition FemtoSHCF.h:241
TH1D * GetNumRe(int el, int em) const
TH1D * GetDenRe(int el, int em) const
TH1D * GetNumIm(int el, int em) const
void Draw(Option_t *opt="")
virtual TString HTMLExtract(Int_t counter=0, TString dir=" ") const
void SetComment(TString comment)
Definition FemtoSHCF.h:141
virtual Long64_t Merge(TCollection *collection)
void Browse(TBrowser *b)
void RecalculateCF(Int_t debugBin=-1, Bool_t suwm=kFALSE)
void SettNormBohr(Double_t bohr)
Definition FemtoSHCF.h:156
void SetNumRe(TH1D **histograms, Bool_t clone=kTRUE)
TH1D * GetRawCFIm(int pos) const
Definition FemtoSHCF.h:247
virtual void FillDenObj(TObject *obj)
void SetScale(Double_t scale)
Definition FemtoSHCF.h:218
TH1D * GetDenIm(int el, int em) const
void SetNormRadius(Double_t radius)
Definition FemtoSHCF.h:151
virtual void Add(const Object *pack)
TH1D * GetCFIm(int el, int em) const
TH3D * GetCovCF() const
void SetNumIm(TH1D **histograms, Bool_t clone=kTRUE)
virtual void FastAdd(const FemtoSHCF *obj)
void AddLabel(TString label)
void Rebin(Int_t ngroup, Option_t *opt="")
void SetNormPurity(Double_t purity)
Definition FemtoSHCF.h:146
Int_t GetLMax() const
Definition FemtoSHCF.h:306
Int_t GetMaxJM() const
Definition FemtoSHCF.h:94
void SetDenRe(TH1D **histograms, Bool_t clone=kTRUE)
virtual TObject * GetSpecial(TString opt) const