BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBunchGaussian.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSBunchGaussian.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSGlobalConstants.hh"
23
24#include "parser/beam.h"
25
26#include "Randomize.hh"
27#include "CLHEP/Matrix/SymMatrix.h"
28#include "CLHEP/Matrix/Vector.h"
29#include "CLHEP/RandomObjects/RandMultiGauss.h"
30#include "CLHEP/Units/PhysicalConstants.h"
31
32#include <cmath>
33#include <vector>
34
36namespace {
38 bool isPositiveDefinite(CLHEP::HepSymMatrix& S)
39 {
40 CLHEP::HepSymMatrix temp (S); // Since diagonalize does not take a const s
41 // we have to copy S.
42 CLHEP::HepMatrix U = diagonalize ( &temp ); // S = U Sdiag U.T()
43 CLHEP::HepSymMatrix D = S.similarityT(U); // D = U.T() S U = Sdiag
44 for (G4int i = 1; i <= S.num_row(); i++)
45 {
46 G4double s2 = D(i,i);
47 if ( s2 <= 0 )
48 {return false;}
49 }
50 return true;
51 }
52}
53
54BDSBunchGaussian::BDSBunchGaussian(const G4String& nameIn):
55 BDSBunch(nameIn),
56 meansGM(CLHEP::HepVector(6)),
57 sigmaGM(CLHEP::HepSymMatrix(6)),
58 gaussMultiGen(nullptr),
59 offsetSampleMean(false),
60 iPartIteration(0)
61{
62 coordinates = {&x0_v, &xp_v, &y0_v, &yp_v, &z0_v, &zp_v, &E_v, &t_v, &weight_v};
63}
64
65BDSBunchGaussian::~BDSBunchGaussian()
66{
67 delete gaussMultiGen;
68}
69
71 const GMAD::Beam& beam,
72 const BDSBunchType& distrType,
73 G4Transform3D beamlineTransformIn,
74 const G4double beamlineSIn)
75{
76 BDSBunch::SetOptions(beamParticle, beam, distrType, beamlineTransformIn, beamlineSIn);
77
78 // require sigmaE to be finite for sigma matrix to be positive definite
79 // set to very small value if 0
80 if (!finiteSigmaE) // set in BDSBunch::SetOptions
81 {sigmaE = 1e-50;}
82 if (!finiteSigmaT)
83 {sigmaT = 1e-50;}
84
85 offsetSampleMean = beam.offsetSampleMean;
87
88 // undo units and redo after the multivariate Gaussian
89 // easier for the Gauss classes where we have sigma^2 in places
90 meansGM[0] = X0 / CLHEP::m;
91 meansGM[1] = Xp0;
92 meansGM[2] = Y0 / CLHEP::m;
93 meansGM[3] = Yp0;
94 meansGM[4] = T0 / CLHEP::s;
95 meansGM[5] = 1;
96}
97
98void BDSBunchGaussian::BeginOfRunAction(G4int numberOfEvents)
99{
100 if (!offsetSampleMean)
101 {return;}
103 for (auto& vec : coordinates)
104 {vec->clear();}
105 iPartIteration = 0;
106
107 PreGenerateEvents(numberOfEvents);
108}
109
110CLHEP::RandMultiGauss* BDSBunchGaussian::CreateMultiGauss(CLHEP::HepRandomEngine& anEngine,
111 const CLHEP::HepVector& mu,
112 CLHEP::HepSymMatrix& sigma)
113{
116 G4cout << __METHOD_NAME__ << "checking 6x6 sigma matrix is positive definite" << G4endl;
117 if (!isPositiveDefinite(sigma))
118 {
119 G4cout << __METHOD_NAME__ << "WARNING bunch generator sigma matrix is not positive definite" << G4endl;
120 G4cout << sigma << G4endl;
121 G4cout << __METHOD_NAME__ << "adding a small error to zero diagonal elements" << G4endl;
122
123 // very small number especially for time, since the sigma goes with the square
124 G4double small_error = 1e-50;
125
126 for (G4int i=0; i<6; i++)
127 {
128 if (sigma[i][i]==0)
129 {sigma[i][i] += small_error;}
130 }
131
132 if (!isPositiveDefinite(sigma))
133 {
134 G4cout << __METHOD_NAME__ << "WARNING bunch generator sigma matrix is still not positive definite" << G4endl;
135 G4cout << sigma << G4endl;
136 G4cout << __METHOD_NAME__ << "adding a small error to all elements" << G4endl;
137 for (G4int i=0; i<6; i++)
138 {
139 for (G4int j=0; j<6; j++)
140 {
141 if (sigma[i][j]==0)
142 {sigma[i][j] += small_error;}
143 }
144 }
145 if (!isPositiveDefinite(sigma))
146 {
147 G4cout << sigma << G4endl;
148 throw BDSException(__METHOD_NAME__, "bunch generator sigma matrix is still not positive definite, giving up");
149 }
150 }
151 }
152
153#ifdef BDSDEBUG
154 G4cout << __METHOD_NAME__ << "mean " << G4endl
155 << mu << G4endl
156 << __METHOD_NAME__ << "sigma " << G4endl
157 << sigma << G4endl;
158#endif
159 G4cout << __METHOD_NAME__ << "confirmed: positive definite matrix" << G4endl;
160 return new CLHEP::RandMultiGauss(anEngine,mu,sigma);
161}
162
164{
165 G4cout << __METHOD_NAME__ << "Pregenerating " << nGenerate << " events." << G4endl;
166 // generate all required primaries first
167 G4double x_a = 0.0, xp_a = 0.0, y_a = 0.0, yp_a = 0.0;
168 G4double z_a = 0.0, zp_a = 0.0, E_a = 0.0, t_a = 0.0;
169
170 for (G4int iParticle = 0; iParticle < nGenerate; ++iParticle)
171 {
173
174 G4double nT = (G4double)iParticle + 1;
175 G4double d = 0;
176 d = fire.x - x_a;
177 x_a = x_a + (d/nT);
178 d = fire.xp - xp_a;
179 xp_a = xp_a + (d/nT);
180 d = fire.y - y_a;
181 y_a = y_a + (d/nT);
182 d = fire.yp - yp_a;
183 yp_a = yp_a + (d/nT);
184 d = fire.z - z_a;
185 z_a = z_a + (d/nT);
186 d = fire.zp - zp_a;
187 zp_a = zp_a + (d/nT);
188 d = fire.totalEnergy - E_a;
189 E_a = E_a + (d/nT);
190 d = fire.T - t_a;
191 t_a = t_a + (d/nT);
192
193 x0_v.push_back(fire.x);
194 xp_v.push_back(fire.xp);
195 y0_v.push_back(fire.y);
196 yp_v.push_back(fire.yp);
197 z0_v.push_back(fire.z);
198 zp_v.push_back(fire.zp);
199 E_v.push_back(fire.totalEnergy);
200 t_v.push_back(fire.T);
201 weight_v.push_back(fire.weight);
202 }
203
204 // Compute difference between sample mean and specified means
205 x_a = x_a - X0;
206 xp_a = xp_a - Xp0;
207 y_a = y_a - Y0;
208 yp_a = yp_a - Yp0;
209 z_a = z_a - Z0;
210 zp_a = zp_a - Zp0;
211 E_a = E_a - E0;
212 t_a = t_a - T0;
213
214 // Offset with different w.r.t. central value
215 for (G4int iParticle = 0; iParticle < nGenerate; ++iParticle)
216 {
217 x0_v[iParticle] -= x_a;
218 xp_v[iParticle] -= xp_a;
219 y0_v[iParticle] -= y_a;
220 yp_v[iParticle] -= yp_a;
221 z0_v[iParticle] -= z_a;
222 zp_v[iParticle] -= zp_a;
223 E_v[iParticle] -= E_a;
224 t_v[iParticle] -= t_a;
225 }
226}
227
229{
231 {
232 // iPartIteration should never exceed the size of each vector.
233 // the units are already correct in the vector of stored coordinates
234 G4double x = x0_v[iPartIteration];
235 G4double xp = xp_v[iPartIteration];
236 G4double y = y0_v[iPartIteration];
237 G4double yp = yp_v[iPartIteration];
238 G4double z = z0_v[iPartIteration];
239 G4double zp = zp_v[iPartIteration];
240 G4double t = t_v[iPartIteration];
241 G4double E = E_v[iPartIteration];
242 G4double weight = weight_v[iPartIteration];
243
245 return BDSParticleCoordsFull(x,y,z,xp,yp,zp,t,S0,E,weight);
246 }
247 else
249}
250
251
253{
254 CLHEP::HepVector v = gaussMultiGen->fire();
255 // unlike other bunch distributions reintroduce units (taken out in set options)
256 G4double x = v[0] * CLHEP::m;
257 G4double xp = v[1];
258 G4double y = v[2] * CLHEP::m;
259 G4double yp = v[3];
260 G4double t = T0;
261 if (finiteSigmaT)
262 {t = v[4];}
263 t *= CLHEP::s;
264 G4double zp = CalculateZp(xp,yp,Zp0);
265 G4double z = Z0;
266 G4double dz = 0;
267 if (finiteSigmaT)
268 {
269 dz = t * CLHEP::c_light;
270 z += dz;
271 }
272 G4double E = E0; // exceptionally left in G4 units
273 if (finiteSigmaE)
274 {E *= v[5];} // only if there's a finite energy spread
275 // cov-matrix is E_bar(1), with sigma fractional - so no units here
276
277 return BDSParticleCoordsFull(x,y,z,xp,yp,zp,t,S0+dz,E,/*weight=*/1.0);
278}
CLHEP::RandMultiGauss * CreateMultiGauss(CLHEP::HepRandomEngine &anEngine, const CLHEP::HepVector &mu, CLHEP::HepSymMatrix &sigma)
std::vector< G4double > z0_v
Holder for pre-calculated coordinates.
std::vector< G4double > y0_v
Holder for pre-calculated coordinates.
std::vector< G4double > xp_v
Holder for pre-calculated coordinates.
virtual BDSParticleCoordsFull GetNextParticleLocal()
std::vector< G4double > t_v
Holder for pre-calculated coordinates.
CLHEP::RandMultiGauss * gaussMultiGen
Randon number generator with sigma matrix and mean.
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
std::vector< std::vector< G4double > * > coordinates
Convenience vector of vectors for clearing up.
std::vector< G4double > weight_v
Holder for pre-calculated coordinates.
virtual void BeginOfRunAction(G4int numberOfEvents)
std::vector< G4double > zp_v
Holder for pre-calculated coordinates.
G4int iPartIteration
Iterator for reading out pre-calculate coordinates.
std::vector< G4double > yp_v
Holder for pre-calculated coordinates.
virtual BDSParticleCoordsFull GetNextParticleLocalCoords()
Fire random number generator and get coordinates. Can be overloaded if required.
std::vector< G4double > E_v
Holder for pre-calculated coordinates.
void PreGenerateEvents(G4int nGenerate)
Pre-generate all the particle coordinates and subtract the sample mean.
G4bool offsetSampleMean
Whether to offset the sample mean.
std::vector< G4double > x0_v
Holder for pre-calculated coordinates.
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
G4double sigmaE
Centre of distributions.
Definition: BDSBunch.hh:173
G4double Yp0
Centre of distributions.
Definition: BDSBunch.hh:166
G4double T0
Centre of distributions.
Definition: BDSBunch.hh:164
G4bool finiteSigmaE
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:189
G4double S0
Centre of distributions.
Definition: BDSBunch.hh:163
G4double sigmaT
Centre of distributions.
Definition: BDSBunch.hh:171
G4double Z0
Centre of distributions.
Definition: BDSBunch.hh:162
G4double X0
Centre of distributions.
Definition: BDSBunch.hh:160
G4double Zp0
Centre of distributions.
Definition: BDSBunch.hh:167
G4double Xp0
Centre of distributions.
Definition: BDSBunch.hh:165
G4double E0
Centre of distributions.
Definition: BDSBunch.hh:168
G4double Y0
Centre of distributions.
Definition: BDSBunch.hh:161
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Definition: BDSBunch.cc:78
G4bool finiteSigmaT
Flags to ignore random number generator in case of no finite E or T.
Definition: BDSBunch.hh:190
static G4double CalculateZp(G4double xp, G4double yp, G4double Zp0)
Calculate zp safely based on other components.
Definition: BDSBunch.cc:317
General exception with possible name of object and message.
Definition: BDSException.hh:35
A set of particle coordinates including energy and weight.
Wrapper for particle definition.
Improve type-safety of native enum data type in C++.
Beam class.
Definition: beam.h:44