Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
CorrFitDumpedPairAna.h
1/*
2 * FemtoNDimMap.h
3 *
4 * Created on: 24 paź 2019
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#ifndef HALFEMTONDIMMAP_H_
10#define HALFEMTONDIMMAP_H_
11
12#include "FemtoFreezoutGenerator.h"
13#include "XMLNode.h"
14
15#include <TObject.h>
16#include <TString.h>
17
18
22class TFile;
23class TChain;
24class TClonesArray;
25namespace Hal {
26 class CorrFitParamsSetup;
27 class CorrFitMapGroupConfig;
28 class Femto1DCF;
29 class Femto3DCF;
30 class FemtoCorrFunc;
31 class FemtoMicroPair;
32 class FemtoPair;
33 class FemtoSHCF;
34 class FemtoWeightGenerator;
35
40 class CorrFitDumpedPairAna : public TObject {
41 std::vector<TString> fUsedBranches;
42 TFile* fFile = {nullptr};
43 TChain* fTree = {nullptr};
44
45 protected:
46 TString fPairFile;
47 Int_t fJobId = {-1};
48 Int_t fMultiplyWeight = {1};
49 Int_t fMultiplyPreprocess = {1};
50 Int_t fMultiplyJobs = {1};
51 Int_t fTotalNumberOfPoints = {0};
52 Int_t fPairThreshold = {0};
53 Int_t fPairsProcessed = {0};
54 Bool_t fIgnoreSing = {kFALSE};
55 Bool_t fImgMom = {kFALSE};
56 FemtoCorrFunc* fTempCF = {nullptr};
57 FemtoPair* fPair = {nullptr};
58 std::vector<FemtoCorrFunc*> fCF;
59 FemtoFreezoutGenerator* fTempGenerator = {nullptr};
60 std::vector<FemtoFreezoutGenerator*> fGenerator;
61 FemtoWeightGenerator* fWeight = {nullptr};
62 CorrFitMapGroupConfig* fGrouping = {nullptr};
63 std::vector<TClonesArray*> fSignalClones;
64 std::vector<TClonesArray*> fBackgroundClones;
65 enum class eDumpCalcMode { kSignalPairs = 0, kSignalBackgroundPairs = 1, kBackgroundPairsOnly = 2 };
66 eDumpCalcMode fMode;
69 Bool_t SaveAsRawArray(TObject* cf, Int_t step);
70 Bool_t ConfigureInput();
71 Bool_t ConfigureRootInput();
72 Bool_t ConfigureListInput();
73 TString FindTreeName(TString name) const;
74 Bool_t ConfigureFromXML();
75 Int_t GetSimStepNo() const { return fMultiplyJobs * fJobId; }
76 virtual void RunSignalPair() = 0;
77 virtual void RunSignalBackgroundPair() = 0;
78 virtual void RunBackgroundPair() = 0;
79 virtual Bool_t IsVertical() const { return kFALSE; }
84 virtual Bool_t ConnectToData() = 0;
85 void LockUnusedBranches();
86 void ConnectToSignal(const std::vector<TString>& branches);
87 void ConnectToBackground(const std::vector<TString>& branches);
92 virtual Bool_t InitCFs() = 0;
96 virtual Bool_t InitGenerators(const std::vector<int>& dims, XMLNode* parameters, const CorrFitParamsSetup& setup) = 0;
101 void SetGenerator(const FemtoFreezoutGenerator& gen) { fTempGenerator = gen.MakeCopy(); };
102
103 public:
104 CorrFitDumpedPairAna(Int_t jobid = -1, Int_t mapsPerAna = -1);
109 void SetCorrFunc(const FemtoCorrFunc& func);
114 void SetWeightGenerator(const FemtoWeightGenerator& weight);
118 void UseBackground() { fMode = eDumpCalcMode::kSignalBackgroundPairs; };
122 void UseBackgroundOnly() { fMode = eDumpCalcMode::kBackgroundPairsOnly; }
126 void UseImgMomenta() { fImgMom = kTRUE; };
127 virtual Bool_t Init();
128 virtual void Run(Int_t maxEvents = -1);
129 virtual void Finish() = 0;
133 void IgnoreSign() { fIgnoreSing = kTRUE; };
138 void SetMultiplyFactorWeight(Int_t m) { fMultiplyWeight = m; };
144 void SetMultiplyFactorPreprocess(Int_t m) { fMultiplyPreprocess = m; };
149 void SetJobID(Int_t jobId) { fJobId = jobId; };
154 void SetPairFile(TString pairFile) { fPairFile = pairFile; };
155 virtual void Print(Option_t* option = "") const;
159 virtual void PreprocessPair() {};
160 virtual void PreprocessMixedPair() {};
161 virtual ~CorrFitDumpedPairAna();
162 ClassDef(CorrFitDumpedPairAna, 1)
163 };
164} // namespace Hal
165#endif /* HALFEMTO_ANALYSIS_FEMTO_CORRFIT_MAPGENERATORS_HALFEMTONDIMMAP_H_ \
166 */
void SetPairFile(TString pairFile)
Bool_t SaveAsRawArray(TObject *cf, Int_t step)
void SetCorrFunc(const FemtoCorrFunc &func)
virtual Bool_t InitCFs()=0
virtual Bool_t InitGenerators(const std::vector< int > &dims, XMLNode *parameters, const CorrFitParamsSetup &setup)=0
void SetWeightGenerator(const FemtoWeightGenerator &weight)
void SetGenerator(const FemtoFreezoutGenerator &gen)
virtual Bool_t ConnectToData()=0
virtual FemtoFreezoutGenerator * MakeCopy() const =0