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,80 @@ R3BCALIFATestGenerator::R3BCALIFATestGenerator()
3334R3BCALIFATestGenerator::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
5740Bool_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!" );
7261
73- for (Int_t i = 0 ; i < fGammaEnergies .size (); i++)
62+ sumBranchingRatios = 0 .;
63+ for (size_t i = 0 ; i < fGammaEnergies .size (); ++i)
7464 {
75- if (fGammaBranchingRatios [i] > 1 )
65+ if (fGammaBranchingRatios [i] > 1 . )
7666 {
7767 R3BLOG (fatal, " gamma branching ratio in position " << i << " larger than 1!" );
7868 }
7969 sumBranchingRatios += fGammaBranchingRatios [i];
8070 }
81- R3BLOG_IF (warn, sumBranchingRatios > 1 , " gamma branching ratio sum larger than 1!" );
71+ R3BLOG_IF (warn, sumBranchingRatios > 1 . , " gamma branching ratio sum larger than 1!" );
8272
83- // Check for particle type
8473 TDatabasePDG* pdgBase = TDatabasePDG::Instance ();
8574 TParticlePDG* particle = pdgBase->GetParticle (fPDGType );
8675 R3BLOG_IF (fatal, !particle, " PDG code " << fPDGType << " not defined." );
8776 fPDGMass = particle->Mass ();
8877 return kTRUE ;
8978}
9079
80+ void R3BCALIFATestGenerator::SampleVertex (Double_t& vx, Double_t& vy, Double_t& vz) const
81+ {
82+ vx = 0 .;
83+ vy = 0 .;
84+ vz = 0 .;
85+
86+ if (fBoxVtxIsSet )
87+ {
88+ vx = gRandom ->Uniform (fX1 , fX2 );
89+ vy = gRandom ->Uniform (fY1 , fY2 );
90+ vz = gRandom ->Uniform (fZ1 , fZ2 );
91+ return ;
92+ }
93+
94+ if (fPointVtxIsSet )
95+ {
96+ if (fDispersionVtxIsSet )
97+ {
98+ vx = gRandom ->Gaus (fX , fDX );
99+ vy = gRandom ->Gaus (fY , fDY );
100+ vz = gRandom ->Uniform (fZ - fDZ , fZ + fDZ );
101+ }
102+ else
103+ {
104+ vx = fX ;
105+ vy = fY ;
106+ vz = fZ ;
107+ }
108+ }
109+ }
110+
91111Bool_t R3BCALIFATestGenerator::ReadEvent (FairPrimaryGenerator* primGen)
92112{
93113 // Generate one event: produce primary particles emitted from one vertex.
@@ -96,11 +116,15 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
96116 // if SetCosTheta() function is used, the distribution will be uniform in
97117 // cos(theta)
98118
99- Double32_t pabs = 0 , phi, pt = 0 , theta = 0 , eta, y, mt, px, py, pz = 0 ;
119+ Double32_t pabs = 0 ., phi = 0 ., pt = 0 ., theta = 0 ., eta = 0 ., y = 0 ., mt = 0 .;
120+ Double32_t px = 0 ., py = 0 ., pz = 0 .;
100121
101122 // Generate particles
102- for (Int_t k = 0 ; k < fMult ; k++ )
123+ for (Int_t k = 0 ; k < fMult ; ++k )
103124 {
125+ Double_t vx = 0 ., vy = 0 ., vz = 0 .;
126+ SampleVertex (vx, vy, vz);
127+
104128 phi = gRandom ->Uniform (fPhiMin , fPhiMax ) * TMath::DegToRad ();
105129
106130 if (fPRangeIsSet )
@@ -118,7 +142,7 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
118142 else if (fEtaRangeIsSet )
119143 {
120144 eta = gRandom ->Uniform (fEtaMin , fEtaMax );
121- theta = 2 * TMath::ATan (TMath::Exp (-eta));
145+ theta = 2 . * TMath::ATan (TMath::Exp (-eta));
122146 }
123147 else if (fYRangeIsSet )
124148 {
@@ -135,7 +159,9 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
135159 pt = pabs * TMath::Sin (theta);
136160 }
137161 else if (fPtRangeIsSet )
162+ {
138163 pz = pt / TMath::Tan (theta);
164+ }
139165 }
140166
141167 px = pt * TMath::Cos (phi);
@@ -151,55 +177,48 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
151177 if (fNuclearDecayChainIsSet )
152178 {
153179 LOG_IF (fatal, fPDGType != 22 ) << " PDG code " << fPDGType << " is not a gamma!" ;
154- for (auto i = 0 ; i < fGammaEnergies .size (); i++ )
180+ for (size_t i = 0 ; i < fGammaEnergies .size (); ++i )
155181 {
156- auto br = gRandom ->Uniform ();
182+ const auto br = gRandom ->Uniform ();
157183 if (br < fGammaBranchingRatios [i])
158184 {
159185 pabs = fGammaEnergies [i];
186+ phi = gRandom ->Uniform (fPhiMin , fPhiMax ) * TMath::DegToRad ();
160187 theta =
161188 acos (gRandom ->Uniform (cos (fThetaMin * TMath::DegToRad ()), cos (fThetaMax * TMath::DegToRad ())));
162189 pz = pabs * TMath::Cos (theta);
163190 pt = pabs * TMath::Sin (theta);
164-
165191 px = pt * TMath::Cos (phi);
166192 py = pt * TMath::Sin (phi);
167193
168- auto gammaMomentum = TMath::Sqrt (px * px + py * py + pz * pz);
194+ const auto gammaMomentum = TMath::Sqrt (px * px + py * py + pz * pz);
169195 px = px * fGammaEnergies [i] / gammaMomentum;
170196 py = py * fGammaEnergies [i] / gammaMomentum;
171197 pz = pz * fGammaEnergies [i] / gammaMomentum;
172- if (sumBranchingRatios <= 1 )
173- break ;
174198
175199 if (fLorentzBoostIsSet )
176200 {
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);
201+ const auto gammaMom = TMath::Sqrt (px * px + py * py + pz * pz);
181202 pz = (pz + fBetaOfEmittingFragment * gammaMom) / fGammaFactor ;
182203 }
183- primGen->AddTrack (fPDGType , px, py, pz, fX , fY , fZ );
204+
205+ primGen->AddTrack (fPDGType , px, py, pz, vx, vy, vz);
206+
207+ if (sumBranchingRatios <= 1 .)
208+ break ;
184209 }
185210 }
186- return kTRUE ;
211+ continue ;
187212 }
188213
189214 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);
215+ {
216+ const Double32_t gammaMomentum = TMath::Sqrt (px * px + py * py + pz * pz);
195217 pz = (pz + fBetaOfEmittingFragment * gammaMomentum) / fGammaFactor ;
196218 }
197219 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 );
220+ {
221+ const Double32_t particleEnergy = TMath::Sqrt (px * px + py * py + pz * pz + fPDGMass * fPDGMass );
203222 pz = (pz + fBetaOfEmittingFragment * particleEnergy) / fGammaFactor ;
204223 }
205224
@@ -208,19 +227,23 @@ Bool_t R3BCALIFATestGenerator::ReadEvent(FairPrimaryGenerator* primGen)
208227 std::ostringstream oss;
209228 oss << " CALIFATestGen: kf=" << fPDGType << " , p=(" << std::fixed << std::setprecision (2 ) << px << " , " << py
210229 << " , " << pz << " ) GeV"
211- << " , x=(" << std::setprecision (1 ) << fX << " , " << fY << " , " << fZ << " ) cm" ;
230+ << " , x=(" << std::setprecision (2 ) << vx << " , " << vy << " , " << vz << " ) cm" ;
231+ R3BLOG (info, oss.str ());
212232 }
213- primGen->AddTrack (fPDGType , px, py, pz, fX , fY , fZ );
233+
234+ primGen->AddTrack (fPDGType , px, py, pz, vx, vy, vz);
214235 }
236+
215237 return kTRUE ;
216238}
217239
218240void R3BCALIFATestGenerator::SetFragmentVelocity (double beta, double dispersion)
219241{
220- // Sets the velocity and gamma factor of the fragment emitting the gammas
242+ // Keep original convention: sample once when the setter is called.
243+ // This matches the current class behavior and avoids unexpected API changes.
221244 R3BLOG (info, " Set a beta of " << beta << " with a sigma dispersion of " << dispersion);
222245 fBetaOfEmittingFragment = gRandom ->Gaus (beta, dispersion);
223- fGammaFactor = TMath::Sqrt (1 - fBetaOfEmittingFragment * fBetaOfEmittingFragment );
246+ fGammaFactor = TMath::Sqrt (1 . - fBetaOfEmittingFragment * fBetaOfEmittingFragment );
224247}
225248
226249void R3BCALIFATestGenerator::SetDecayChainPoint (double gammaEnergy, double branchingRatio)
@@ -229,4 +252,5 @@ void R3BCALIFATestGenerator::SetDecayChainPoint(double gammaEnergy, double branc
229252 fGammaEnergies .push_back (gammaEnergy / 1000 .); // GeV
230253 fGammaBranchingRatios .push_back (branchingRatio);
231254}
255+
232256ClassImp (R3BCALIFATestGenerator)
0 commit comments