Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
CorrFitPairGeneratorSimple.cxx
1/*
2 * CorrFitPairGeneratorSimple.cxx
3 *
4 * Created on: 27 lut 2024
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "CorrFitPairGeneratorSimple.h"
11#include "CorrFitMapGroupConfig.h"
12
13#include "Array.h"
14#include "Cout.h"
15#include "FemtoMiniPair.h"
16
17#include <RtypesCore.h>
18#include <TClonesArray.h>
19#include <TLorentzVector.h>
20#include <TObjArray.h>
21#include <TRandom.h>
22
23namespace Hal {
24
25 void CorrFitPairGeneratorSimple::GenerateEvent() {
26 if (fLimitsN[0] > fBinLimit) return;
27 fLimitsN[0] += 100;
28 Double_t min = fGrouping.GetMin();
29 Double_t max = fGrouping.GetMax();
30 Double_t scale = 1.0;
31 if (fFrame == Femto::EKinematics::kLCMS) { scale = 0.5; }
32 min = min * scale;
33 max = max * scale;
34 TLorentzVector p1, p2;
35 Double_t x, y, z;
36 FemtoMicroPair* pair = nullptr;
37 auto build = [&]() {
38 p1.SetXYZM(x, y, z, fM1);
39 p2.SetXYZM(-x, -y, -z, fM2);
40 p1.Boost(0.2, 0, 0.0);
41 p2.Boost(0.2, 0, 0.0);
42 pair->SetMomenta1(p1.X(), p1.Y(), p1.Z(), p1.E());
43 pair->SetTrueMomenta1(p1.X(), p1.Y(), p1.Z(), p1.E());
44 pair->SetMomenta2(p2.X(), p2.Y(), p2.Z(), p2.E());
45 pair->SetTrueMomenta2(p2.X(), p2.Y(), p2.Z(), p2.E());
46 pair->SetPdg1(fPid1);
47 pair->SetPdg2(fPid2);
48 };
49
50
51 if (fGrouping.GroupByKStar()) {
52 const Int_t defSize = 100;
53 for (int bin = 0; bin < fGrouping.GetNbins(); bin++) {
54 const Double_t kstar = fCentersX[bin] * scale;
55 fSignalPairs[bin]->ExpandCreateFast(defSize);
56 for (int i = 0; i < defSize; i++) {
57 pair = (FemtoMicroPair*) fSignalPairs[bin]->UncheckedAt(i);
58 gRandom->Sphere(x, y, z, kstar);
59 build();
60 }
61 }
62 } else {
63 const Int_t defSize = 1;
64 const int binY = fCentersY.GetSize();
65 const int binX = fCentersX.GetSize();
66 for (int bin = 0; bin < fGrouping.GetNbins(); bin++) {
67 fSignalPairs[bin]->ExpandCreateFast(defSize * binY * binX);
68 int count = 0;
69 for (int biny = 0; biny < binY; biny++) {
70 for (int binx = 0; binx < binX; binx++) {
71 pair = (FemtoMicroPair*) fSignalPairs[bin]->UncheckedAt(count++);
72 z = fCentersZ[bin] * scale;
73 y = fCentersY[biny] * scale;
74 x = fCentersX[binx] * scale;
75 if (gRandom->Uniform() < 0.5) x = -x;
76 if (gRandom->Uniform() < 0.5) y = -y;
77 if (gRandom->Uniform() < 0.5) z = -z;
78 build();
79 }
80 }
81 }
82 } /* finished 3d generating*/
83 }
84
85 Bool_t CorrFitPairGeneratorSimple::Init() {
86 auto val = CorrFitPairGenerator::Init();
87 if (fFrame != Hal::Femto::EKinematics::kPRF)
88 Hal::Cout::PrintInfo(Form("CorrFitPairGeneratorSimple works correctly only with PRF frame"), EInfo::kError);
89 return val;
90 }
91
92} /* namespace Hal */
Int_t GetSize() const
Definition Array.h:50
static void PrintInfo(TString text, Hal::EInfo status)
Definition Cout.cxx:370