BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSPhysicsLaserWire.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 "BDSLaserCompton.hh"
20#include "BDSPhysicsLaserWire.hh"
21
22#include "globals.hh" // geant4 types / globals
23#include "G4Electron.hh"
24#include "G4Gamma.hh"
25#include "G4OpticalPhoton.hh"
26#include "G4ParticleDefinition.hh"
27#include "G4Positron.hh"
28#include "G4ProcessManager.hh"
29#include "G4Version.hh"
30
31BDSPhysicsLaserWire::BDSPhysicsLaserWire():
32 G4VPhysicsConstructor("BDSPhysicsLaserWire")
33{;}
34
35BDSPhysicsLaserWire::~BDSPhysicsLaserWire()
36{;}
37
39{
40 G4Electron::ElectronDefinition();
41 G4Positron::PositronDefinition();
42 G4Gamma::Gamma();
43 G4OpticalPhoton::OpticalPhotonDefinition();
44}
45
47{
48 if (Activated())
49 {return;}
50
51#if G4VERSION_NUMBER > 1029
52 auto aParticleIterator = GetParticleIterator();
53#endif
54 aParticleIterator->reset();
55
56 BDSLaserCompton* lwProcess = new BDSLaserCompton();
57
58 while((*aParticleIterator)())
59 {
60 G4ParticleDefinition* particle = aParticleIterator->value();
61 G4ProcessManager* pmanager = particle->GetProcessManager();
62 G4String particleName = particle->GetParticleName();
63
64 if (particleName == "e-")
65 {
66 pmanager->AddProcess(lwProcess);
67 pmanager->SetProcessOrderingToLast(lwProcess,idxPostStep);
68 }
69
70 if (particleName == "e+")
71 {
72 pmanager->AddProcess(lwProcess);
73 pmanager->SetProcessOrderingToLast(lwProcess,idxPostStep);
74 }
75 }
76
78}
Laser compton scattering process to achieve a laserwire.
virtual void ConstructParticle()
Construct electrons and positrons and the photon.
virtual void ConstructProcess()
Construct and attache BDSIM laserwire process.
G4bool Activated() const
Get whether this instance has been activated.
Definition: BDSSingleUse.hh:37
void SetActivated()
Flag this instance as activated for later querying.
Definition: BDSSingleUse.hh:40