BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSLinkPrimaryGeneratorAction.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 "BDSBunch.hh"
20#include "BDSBunchSixTrackLink.hh"
21#include "BDSDebug.hh"
22#include "BDSEventInfo.hh"
23#include "BDSException.hh"
24#include "BDSExtent.hh"
25#include "BDSIonDefinition.hh"
26#include "BDSLinkDetectorConstruction.hh"
27#include "BDSLinkEventInfo.hh"
28#include "BDSLinkPrimaryGeneratorAction.hh"
29#include "BDSLinkRegistry.hh"
30#include "BDSPrimaryVertexInformation.hh"
31#include "BDSRandom.hh"
32
33#include "G4Event.hh"
34#include "G4IonTable.hh"
35#include "G4ParticleGun.hh"
36#include "G4Types.hh"
37
39 int* currentElementIndexIn,
40 BDSLinkDetectorConstruction* constructionIn,
41 G4bool debugIn):
42 bunch(bunchIn),
43 currentElementIndex(currentElementIndexIn),
44 construction(constructionIn),
45 debug(debugIn),
46 particleGun(nullptr)
47{
48 particleGun = new G4ParticleGun(1); // 1-particle gun
49
50 particleGun->SetParticleMomentumDirection(G4ThreeVector(0.,0.,1.));
51 particleGun->SetParticlePosition(G4ThreeVector());
52 particleGun->SetParticleTime(0);
53}
54
55BDSLinkPrimaryGeneratorAction::~BDSLinkPrimaryGeneratorAction()
56{
57 delete particleGun;
58}
59
61{
62 // always save seed state in output
63 BDSLinkEventInfo* eventInfo = new BDSLinkEventInfo();
64 anEvent->SetUserInformation(eventInfo);
65 eventInfo->SetSeedStateAtStart(BDSRandom::GetSeedState());
66
68 try
69 {
70 coords = bunch->GetNextParticleLocal();
71 auto bunchSTL = dynamic_cast<BDSBunchSixTrackLink*>(bunch);
72 if (bunchSTL)
73 {
74 eventInfo->externalParticleIDofPrimary = bunchSTL->CurrentExternalParticleID();
75 eventInfo->externalParentIDofPrimary = bunchSTL->CurrentExternalParentID();
76 }
77 }
78 catch (const BDSException& exception)
79 {// we couldn't safely generate a particle -> abort
80 // could be because of user input file
81 anEvent->SetEventAborted();
82 G4cout << exception.what() << G4endl;
83 G4cout << "Aborting this event (#" << anEvent->GetEventID() << ")" << G4endl;
84 return;
85 }
86
88 auto lr = construction->LinkRegistry();
89 const G4Transform3D tr = lr->Transform(*currentElementIndex);
90 if (lr->NoRotation(*currentElementIndex))
91 {
92 if (debug)
93 {
94 G4cout << "PGA: Coords before " << coords;
95 G4cout << "Offset " << tr.getTranslation() << G4endl;
96 }
97 BDSParticleCoords cgf = coords.ApplyOffset(tr.getTranslation());
98 cg = BDSParticleCoordsFullGlobal(coords, cgf);
99 if (debug)
100 {G4cout << "Coords after " << cg.global;}
101 }
102 else
103 {
105 }
106
107 particleGun->SetParticleDefinition(bunch->ParticleDefinition()->ParticleDefinition());
108
109 // always update the charge - ok for normal particles; fixes purposively specified ions.
110 particleGun->SetParticleCharge(bunch->ParticleDefinition()->Charge());
111
112 // check that kinetic energy is positive and finite anyway and abort if not.
113 // get the mass from the beamParticle as this takes into account any electrons
114 G4double EK = cg.local.totalEnergy - bunch->ParticleDefinition()->Mass();
115 if (EK <= 0)
116 {
117 G4cout << __METHOD_NAME__ << "Event #" << anEvent->GetEventID()
118 << " - Particle kinetic energy smaller than 0! "
119 << "This will not be tracked." << G4endl;
120 anEvent->SetEventAborted();
121 return;
122 }
123
124 // check the coordinates are valid
125 if (!worldExtent.Encompasses(cg.global))
126 {
127 G4cerr << __METHOD_NAME__ << "point: " << cg.global
128 << "mm lies outside the world volume with extent ("
129 << worldExtent << " - event aborted!" << G4endl << G4endl;
130 anEvent->SetEventAborted();
131 }
132
133#ifdef BDSDEBUG
134 G4cout << __METHOD_NAME__ << coords << G4endl;
135#endif
136
137 G4ThreeVector PartMomDir(cg.global.xp,cg.global.yp,cg.global.zp);
138 G4ThreeVector PartPosition(cg.global.x,cg.global.y,cg.global.z);
139
140 particleGun->SetParticlePosition(PartPosition);
141 particleGun->SetParticleEnergy(EK);
142 particleGun->SetParticleMomentumDirection(PartMomDir);
143 particleGun->SetParticleTime(cg.global.T);
144
145 particleGun->GeneratePrimaryVertex(anEvent);
146
147 // set the weight
148 auto vertex = anEvent->GetPrimaryVertex();
149 vertex->SetWeight(cg.local.weight);
150 //vertex->Print();
151
152 // associate full set of coordinates with vertex for writing to output after event
153 //vertex->SetUserInformation(new BDSPrimaryVertexInformation(coords,
154 // bunch->ParticleDefinition()));
155
156#ifdef BDSDEBUG
157 vertex->Print();
158#endif
159}
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
virtual const BDSParticleDefinition * ParticleDefinition() const
Access the beam particle definition.
Definition: BDSBunch.hh:95
virtual BDSParticleCoordsFull GetNextParticleLocal()
Definition: BDSBunch.cc:244
void SetSeedStateAtStart(const G4String &seedStateAtStartIn)
Setters.
Definition: BDSEventInfo.hh:54
General exception with possible name of object and message.
Definition: BDSException.hh:35
const char * what() const noexcept override
Override message in std::exception.
Definition: BDSException.hh:55
G4bool Encompasses(const G4ThreeVector &point) const
Return whether the extent encompasses the point. True if point lies inside the extent.
Definition: BDSExtent.cc:190
Construction of the geometry in the case of a link model.
BDSLinkRegistry * LinkRegistry() const
Accessor.
Simple extension to cache extra variables through an event.
BDSBunch * bunch
BDSIM particle generator.
virtual void GeneratePrimaries(G4Event *)
Main interface for Geant4. Prepare primary(ies) for the event.
int * currentElementIndex
External integer for which element to track in.
G4ParticleGun * particleGun
Geant4 particle gun that creates single particles.
BDSExtent worldExtent
World extent that particle coordinates are checked against to ensure they're inside it.
BDSLinkPrimaryGeneratorAction(BDSBunch *bunchIn, int *currentElementIndexIn, BDSLinkDetectorConstruction *constructionIn, G4bool debugIn=false)
Bunch must have a valid particle definition (ie not nullptr).
BDSLinkDetectorConstruction * construction
Cache of detector construction for link registry of transforms.
A set of particle coordinates in both local and global.
A set of particle coordinates including energy and weight.
A set of particle coordinates.
BDSParticleCoords ApplyTransform(const G4Transform3D &transform) const
Apply a transform to the coordinates and return a copy of them transformed.
BDSParticleCoords ApplyOffset(const G4ThreeVector &offset) const
Apply an offset to the spatial coordinates only and return a copy.
G4double Mass() const
Accessor.
G4double Charge() const
Accessor.
G4ParticleDefinition * ParticleDefinition() const
Accessor.