BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBOptrMultiParticleChangeCrossSection.hh
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#ifndef BDSBOPTRMULTIPARTICLECHANGECROSSSECTION_H
20#define BDSBOPTRMULTIPARTICLECHANGECROSSSECTION_H
21
22#include "G4Version.hh"
23#if G4VERSION_NUMBER > 1009
24
25#include "G4VBiasingOperator.hh"
27class G4ParticleDefinition;
28
29#include <map>
30
42class BDSBOptrMultiParticleChangeCrossSection: public G4VBiasingOperator
43{
44public:
47
48 void AddParticle(const G4String& particleName);
49 void SetBias(const G4String& biasObjectName,
50 const G4String& particleName,
51 const G4String& process,
52 G4double dBias,
53 G4int iPrimary);
54 void StartTracking(const G4Track* track) override;
55
56private:
57 virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation(const G4Track* track,
58 const G4BiasingProcessInterface* callingProcess) override;
59 // -- Methods not used:
60 virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation(const G4Track*, const G4BiasingProcessInterface*) override
61 {return 0;}
62 virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation(const G4Track*, const G4BiasingProcessInterface*) override
63 {return 0;}
64 virtual void OperationApplied(const G4BiasingProcessInterface* callingProcess,
65 G4BiasingAppliedCase biasingCase,
66 G4VBiasingOperation* occurenceOperationApplied,
67 G4double weightForOccurenceInteraction,
68 G4VBiasingOperation* finalStateOperationApplied,
69 const G4VParticleChange* particleChangeProduced) override;
70 // prevent compiler warning (since second G4VBiasingOperator::OperationApplied is hidden)
71 using G4VBiasingOperator::OperationApplied;
72
73 std::map<const G4ParticleDefinition*, BDSBOptrChangeCrossSection*> fBOptrForParticle;
74 std::vector<const G4ParticleDefinition*> fParticlesToBias;
75 BDSBOptrChangeCrossSection* fCurrentOperator;
76
78 G4bool debug;
79};
80
81#endif
82
83#endif
Operator that changes cross section of a process for a particle.
G4int fnInteractions
Count number of biased interations for current track.