Skip to content

Commit 3f9ed9a

Browse files
Mrunmoy-Jenajose-luis-rs
authored andcommitted
feat(r3bgen): Added INCL-like vertex smearing in CALIFATestGenerator
clang-formatting applied clang-formatting fixed for R3BCALIFATestGenerator.cxx Fixed headers back to original Added angular spread of the outgoing fragment and sample theta from it for accurate Doppler correction aaply clang-format
1 parent 31e256e commit 3f9ed9a

2 files changed

Lines changed: 272 additions & 121 deletions

File tree

r3bgen/R3BCALIFATestGenerator.cxx

Lines changed: 138 additions & 60 deletions
Original file line numberDiff line numberDiff line change
@@ -20,6 +20,7 @@
2020
#include <TMath.h>
2121
#include <TParticlePDG.h>
2222
#include <TRandom.h>
23+
2324
#include <iomanip>
2425
#include <iostream>
2526
#include <sstream>
@@ -33,61 +34,131 @@ R3BCALIFATestGenerator::R3BCALIFATestGenerator()
3334
R3BCALIFATestGenerator::R3BCALIFATestGenerator(Int_t pdgid, Int_t mult)
3435
: fPDGType(pdgid)
3536
, fMult(mult)
36-
, fPDGMass(0)
37-
, fPtMin(0)
38-
, fPtMax(0)
39-
, fEtaMin(0)
40-
, fEtaMax(0)
41-
, fYMin(0)
42-
, fYMax(0)
43-
, fPMin(0)
44-
, fPMax(0)
45-
, fX(0)
46-
, fY(0)
47-
, fZ(0)
48-
, fX1(0)
49-
, fY1(0)
50-
, fZ1(0)
51-
, fX2(0)
52-
, fY2(0)
53-
, fZ2(0)
5437
{
5538
}
5639

5740
Bool_t R3BCALIFATestGenerator::Init()
5841
{
5942
// Initialize generator
60-
R3BLOG_IF(fatal, fPhiMax - fPhiMin > 360, "phi range is too wide: " << fPhiMin << "<phi<" << fPhiMax);
43+
R3BLOG_IF(fatal, fPhiMax - fPhiMin > 360., "phi range is too wide: " << fPhiMin << "<phi<" << fPhiMax);
6144
R3BLOG_IF(fatal, fPRangeIsSet && fPtRangeIsSet, "Cannot set P and Pt ranges simultaneously");
6245
R3BLOG_IF(fatal, fPRangeIsSet && fYRangeIsSet, "Cannot set P and Y ranges simultaneously");
6346

6447
if ((fThetaRangeIsSet && fYRangeIsSet) || (fThetaRangeIsSet && fEtaRangeIsSet) || (fYRangeIsSet && fEtaRangeIsSet))
6548
{
6649
R3BLOG(fatal, "Cannot set Y, Theta or Eta ranges simultaneously");
6750
}
68-
R3BLOG_IF(fatal, fPointVtxIsSet && fBoxVtxIsSet, "Cannot set point and box vertices simultaneously");
51+
52+
R3BLOG_IF(fatal,
53+
fBoxVtxIsSet && (fPointVtxIsSet || fDispersionVtxIsSet),
54+
"Cannot set explicit box vertex together with SetXYZ/SetDxDyDz vertex conventions");
55+
R3BLOG_IF(fatal,
56+
fDispersionVtxIsSet && !fPointVtxIsSet,
57+
"SetDxDyDz requires a central vertex to be defined first via SetXYZ");
6958

7059
// CALIFA specifics
71-
R3BLOG_IF(fatal, fBetaOfEmittingFragment > 1, "beta of fragment larger than 1!");
60+
R3BLOG_IF(fatal, fBetaOfEmittingFragment > 1., "beta of fragment larger than 1!");
61+
R3BLOG_IF(fatal,
62+
fBeamAngularDivIsSet && !fLorentzBoostIsSet,
63+
"SetBeamAngularDivergence has no effect without SetLorentzBoost");
7264

73-
for (Int_t i = 0; i < fGammaEnergies.size(); i++)
65+
sumBranchingRatios = 0.;
66+
for (size_t i = 0; i < fGammaEnergies.size(); ++i)
7467
{
75-
if (fGammaBranchingRatios[i] > 1)
68+
if (fGammaBranchingRatios[i] > 1.)
7669
{
7770
R3BLOG(fatal, "gamma branching ratio in position " << i << " larger than 1!");
7871
}
7972
sumBranchingRatios += fGammaBranchingRatios[i];
8073
}
81-
R3BLOG_IF(warn, sumBranchingRatios > 1, "gamma branching ratio sum larger than 1!");
74+
R3BLOG_IF(warn, sumBranchingRatios > 1., "gamma branching ratio sum larger than 1!");
8275

83-
// Check for particle type
8476
TDatabasePDG* pdgBase = TDatabasePDG::Instance();
8577
TParticlePDG* particle = pdgBase->GetParticle(fPDGType);
8678
R3BLOG_IF(fatal, !particle, "PDG code " << fPDGType << " not defined.");
8779
fPDGMass = particle->Mass();
8880
return kTRUE;
8981
}
9082

83+
void R3BCALIFATestGenerator::SampleVertex(Double_t& vx, Double_t& vy, Double_t& vz) const
84+
{
85+
vx = 0.;
86+
vy = 0.;
87+
vz = 0.;
88+
89+
if (fBoxVtxIsSet)
90+
{
91+
vx = gRandom->Uniform(fX1, fX2);
92+
vy = gRandom->Uniform(fY1, fY2);
93+
vz = gRandom->Uniform(fZ1, fZ2);
94+
return;
95+
}
96+
97+
if (fPointVtxIsSet)
98+
{
99+
if (fDispersionVtxIsSet)
100+
{
101+
vx = gRandom->Gaus(fX, fDX);
102+
vy = gRandom->Gaus(fY, fDY);
103+
vz = gRandom->Uniform(fZ - fDZ, fZ + fDZ);
104+
}
105+
else
106+
{
107+
vx = fX;
108+
vy = fY;
109+
vz = fZ;
110+
}
111+
}
112+
}
113+
114+
void R3BCALIFATestGenerator::SampleBeamDirection(Double_t& nx, Double_t& ny, Double_t& nz) const
115+
{
116+
nx = 0.;
117+
ny = 0.;
118+
nz = 1.;
119+
120+
if (!fBeamAngularDivIsSet)
121+
return;
122+
123+
// Sample projected angles in the x-z and y-z planes from a Gaussian
124+
// centered on zero (nominal beam axis)
125+
const Double_t thetaX = gRandom->Gaus(0., fBeamDivX);
126+
const Double_t thetaY = gRandom->Gaus(0., fBeamDivY);
127+
nx = TMath::Sin(thetaX);
128+
ny = TMath::Sin(thetaY);
129+
// Ensure unit vector; clamp argument to avoid sqrt of negative from
130+
// floating-point rounding when divergence angles are very large.
131+
const Double_t nz2 = 1. - nx * nx - ny * ny;
132+
nz = (nz2 > 0.) ? TMath::Sqrt(nz2) : 0.;
133+
}
134+
135+
void R3BCALIFATestGenerator::ApplyLorentzBoost(Double_t& px,
136+
Double_t& py,
137+
Double_t& pz,
138+
Double_t nx,
139+
Double_t ny,
140+
Double_t nz) const
141+
{
142+
// General Lorentz boost along unit vector n = (nx, ny, nz):
143+
// E' = gamma * (E + beta * (n . p))
144+
// p' = p + n * [(gamma - 1) * (n . p) + gamma * beta * E]
145+
//
146+
// fGammaFactor stores sqrt(1 - beta^2) = 1/gamma, matching the
147+
// convention used in SetFragmentVelocity.
148+
149+
const Double_t E = (fPDGType == 22) ? TMath::Sqrt(px * px + py * py + pz * pz)
150+
: TMath::Sqrt(px * px + py * py + pz * pz + fPDGMass * fPDGMass);
151+
152+
const Double_t gamma = 1. / fGammaFactor;
153+
const Double_t beta = fBetaOfEmittingFragment;
154+
const Double_t pDotN = px * nx + py * ny + pz * nz;
155+
const Double_t coeff = (gamma - 1.) * pDotN + gamma * beta * E;
156+
157+
px += nx * coeff;
158+
py += ny * coeff;
159+
pz += nz * coeff;
160+
}
161+
91162
Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
92163
{
93164
// Generate one event: produce primary particles emitted from one vertex.
@@ -96,11 +167,21 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
96167
// if SetCosTheta() function is used, the distribution will be uniform in
97168
// cos(theta)
98169

99-
Double32_t pabs = 0, phi, pt = 0, theta = 0, eta, y, mt, px, py, pz = 0;
170+
Double32_t pabs = 0., phi = 0., pt = 0., theta = 0., eta = 0., y = 0., mt = 0.;
171+
Double32_t px = 0., py = 0., pz = 0.;
172+
173+
// Sample the boost direction once for the whole event: all de-excitation
174+
// photons from a single fragment share the same beam direction.
175+
Double_t beamNx = 0., beamNy = 0., beamNz = 1.;
176+
if (fLorentzBoostIsSet)
177+
SampleBeamDirection(beamNx, beamNy, beamNz);
100178

101179
// Generate particles
102-
for (Int_t k = 0; k < fMult; k++)
180+
for (Int_t k = 0; k < fMult; ++k)
103181
{
182+
Double_t vx = 0., vy = 0., vz = 0.;
183+
SampleVertex(vx, vy, vz);
184+
104185
phi = gRandom->Uniform(fPhiMin, fPhiMax) * TMath::DegToRad();
105186

106187
if (fPRangeIsSet)
@@ -118,7 +199,7 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
118199
else if (fEtaRangeIsSet)
119200
{
120201
eta = gRandom->Uniform(fEtaMin, fEtaMax);
121-
theta = 2 * TMath::ATan(TMath::Exp(-eta));
202+
theta = 2. * TMath::ATan(TMath::Exp(-eta));
122203
}
123204
else if (fYRangeIsSet)
124205
{
@@ -135,7 +216,9 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
135216
pt = pabs * TMath::Sin(theta);
136217
}
137218
else if (fPtRangeIsSet)
219+
{
138220
pz = pt / TMath::Tan(theta);
221+
}
139222
}
140223

141224
px = pt * TMath::Cos(phi);
@@ -151,76 +234,70 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
151234
if (fNuclearDecayChainIsSet)
152235
{
153236
LOG_IF(fatal, fPDGType != 22) << "PDG code " << fPDGType << " is not a gamma!";
154-
for (auto i = 0; i < fGammaEnergies.size(); i++)
237+
for (size_t i = 0; i < fGammaEnergies.size(); ++i)
155238
{
156-
auto br = gRandom->Uniform();
239+
const auto br = gRandom->Uniform();
157240
if (br < fGammaBranchingRatios[i])
158241
{
159242
pabs = fGammaEnergies[i];
243+
phi = gRandom->Uniform(fPhiMin, fPhiMax) * TMath::DegToRad();
160244
theta =
161245
acos(gRandom->Uniform(cos(fThetaMin * TMath::DegToRad()), cos(fThetaMax * TMath::DegToRad())));
162246
pz = pabs * TMath::Cos(theta);
163247
pt = pabs * TMath::Sin(theta);
164-
165248
px = pt * TMath::Cos(phi);
166249
py = pt * TMath::Sin(phi);
167250

168-
auto gammaMomentum = TMath::Sqrt(px * px + py * py + pz * pz);
251+
const auto gammaMomentum = TMath::Sqrt(px * px + py * py + pz * pz);
169252
px = px * fGammaEnergies[i] / gammaMomentum;
170253
py = py * fGammaEnergies[i] / gammaMomentum;
171254
pz = pz * fGammaEnergies[i] / gammaMomentum;
172-
if (sumBranchingRatios <= 1)
173-
break;
174255

175256
if (fLorentzBoostIsSet)
176257
{
177-
// Lorentz transformation Pz(lab) = gamma * Pz(cm) + gamma * beta * E
178-
// As each Lorentz transformation can be performed sequentially,
179-
// we can separate the gamma factor corresponding to each direction
180-
auto gammaMom = TMath::Sqrt(px * px + py * py + pz * pz);
181-
pz = (pz + fBetaOfEmittingFragment * gammaMom) / fGammaFactor;
258+
ApplyLorentzBoost(px, py, pz, beamNx, beamNy, beamNz);
182259
}
183-
primGen->AddTrack(fPDGType, px, py, pz, fX, fY, fZ);
260+
261+
primGen->AddTrack(fPDGType, px, py, pz, vx, vy, vz);
262+
263+
if (sumBranchingRatios <= 1.)
264+
break;
184265
}
185266
}
186-
return kTRUE;
267+
continue;
187268
}
188269

189-
if (fPDGType == 22 && fLorentzBoostIsSet)
190-
{ // for gamma-rays
191-
// Lorentz transformation Pz(lab) = gamma * Pz(cm) + gamma * beta * E
192-
// As each Lorentz transformation can be performed sequentially,
193-
// we can separate the gamma factor corresponding to each direction
194-
Double32_t gammaMomentum = TMath::Sqrt(px * px + py * py + pz * pz);
195-
pz = (pz + fBetaOfEmittingFragment * gammaMomentum) / fGammaFactor;
196-
}
197-
else if (fLorentzBoostIsSet)
198-
{ // for any massive particle
199-
// Lorentz transformation Pz(lab) = gamma * Pz(cm) + gamma * beta * E
200-
// As each Lorentz transformation can be performed sequentially,
201-
// we can separate the gamma factor corresponding to each direction
202-
Double32_t particleEnergy = TMath::Sqrt(px * px + py * py + pz * pz + fPDGMass * fPDGMass);
203-
pz = (pz + fBetaOfEmittingFragment * particleEnergy) / fGammaFactor;
270+
if (fLorentzBoostIsSet)
271+
{
272+
ApplyLorentzBoost(px, py, pz, beamNx, beamNy, beamNz);
204273
}
205274

206275
if (fDebug)
207276
{
208277
std::ostringstream oss;
209278
oss << "CALIFATestGen: kf=" << fPDGType << ", p=(" << std::fixed << std::setprecision(2) << px << ", " << py
210279
<< ", " << pz << ") GeV"
211-
<< ", x=(" << std::setprecision(1) << fX << ", " << fY << ", " << fZ << ") cm";
280+
<< ", x=(" << std::setprecision(2) << vx << ", " << vy << ", " << vz << ") cm";
281+
if (fLorentzBoostIsSet && fBeamAngularDivIsSet)
282+
{
283+
oss << ", beam_n=(" << std::setprecision(4) << beamNx << ", " << beamNy << ", " << beamNz << ")";
284+
}
285+
R3BLOG(info, oss.str());
212286
}
213-
primGen->AddTrack(fPDGType, px, py, pz, fX, fY, fZ);
287+
288+
primGen->AddTrack(fPDGType, px, py, pz, vx, vy, vz);
214289
}
290+
215291
return kTRUE;
216292
}
217293

218294
void R3BCALIFATestGenerator::SetFragmentVelocity(double beta, double dispersion)
219295
{
220-
// Sets the velocity and gamma factor of the fragment emitting the gammas
296+
// Keep original convention: sample once when the setter is called.
297+
// This matches the current class behavior and avoids unexpected API changes.
221298
R3BLOG(info, "Set a beta of " << beta << " with a sigma dispersion of " << dispersion);
222299
fBetaOfEmittingFragment = gRandom->Gaus(beta, dispersion);
223-
fGammaFactor = TMath::Sqrt(1 - fBetaOfEmittingFragment * fBetaOfEmittingFragment);
300+
fGammaFactor = TMath::Sqrt(1. - fBetaOfEmittingFragment * fBetaOfEmittingFragment);
224301
}
225302

226303
void R3BCALIFATestGenerator::SetDecayChainPoint(double gammaEnergy, double branchingRatio)
@@ -229,4 +306,5 @@ void R3BCALIFATestGenerator::SetDecayChainPoint(double gammaEnergy, double branc
229306
fGammaEnergies.push_back(gammaEnergy / 1000.); // GeV
230307
fGammaBranchingRatios.push_back(branchingRatio);
231308
}
309+
232310
ClassImp(R3BCALIFATestGenerator)

0 commit comments

Comments
 (0)