Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAScavengerProcess.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27#ifndef G4DNASCAVENGERPROCESS_HH
28#define G4DNASCAVENGERPROCESS_HH
29
30#include "G4VITProcess.hh"
31#include "G4MoleculeTable.hh"
32
38{
39 public:
42 explicit G4DNAScavengerProcess(const G4String& aName,
43 const G4DNABoundingBox& box,
45 ~G4DNAScavengerProcess() override;
48 void StartTracking(G4Track*) override;
49 void SetReaction(MolType, Data* pData);
50
51 public:
52 void BuildPhysicsTable(const G4ParticleDefinition&) override;
53
55 const G4Track& track, G4double previousStepSize,
56 G4ForceCondition* condition) override;
57
58 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
59
61 G4ForceCondition*) override
62 {
63 return -1.0;
64 }
65
66 G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override
67 {
68 return nullptr;
69 }
70
71 // no operation in AlongStepDoIt
74 G4GPILSelection*) override
75 {
76 return -1.0;
77 }
78
79 // no operation in AlongStepDoIt
81 {
82 return nullptr;
83 }
84
85 protected:
93
94 protected:
99 std::map<MolType /*MolConf*/, std::map<MolType /*molConfMat*/, Data*>>
101 std::vector<MolType> fpMaterialVector;
108};
109#endif // FLASH1_G4DNASCAVENGERPROCESS_HH
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4ProcessType
@ fUserDefined
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
const G4MolecularConfiguration * MolType
void StartTracking(G4Track *) override
std::map< MolType, std::map< MolType, Data * > > fConfMap
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
const G4DNABoundingBox * fpBoundingBox
const G4MolecularConfiguration * fpMolecularConfiguration
G4DNAScavengerProcess(const G4String &aName, const G4DNABoundingBox &box, G4ProcessType type=fUserDefined)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
std::vector< MolType > fpMaterialVector
G4DNAScavengerMaterial * fpScavengerMaterial
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4DNAScavengerProcess(const G4DNAScavengerProcess &)=delete
G4DNAScavengerProcess & operator=(const G4DNAScavengerProcess &)=delete
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
void SetReaction(MolType, Data *pData)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
const G4DNAMolecularReactionData Data
G4MolecularConfiguration * GetConfiguration(const G4String &, bool mustExist=true)
static G4MoleculeTable * Instance()
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)