Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPolyNucleotideReactionProcess.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/// file: G4DNAPolyNucleotideReactionProcess.hh
27/// brief: This file handls reaction process with DNA geometry.
28#ifndef G4DNAPolyNucleotideReactionProcess_hh
29#define G4DNAPolyNucleotideReactionProcess_hh
32#include <variant>
34class G4DNAComponentNode;
35class G4VDNAHitModel;
36
38{
39 public:
41 const G4String& aName = "DNAStaticMoleculeReactionProcess",
42 G4int verbosityLevel = 0);
44
45 inline void SetDNADamageReactionModel(G4VDNAHitModel* pModel);
46
47 G4bool IsApplicable(const G4ParticleDefinition&) override { return true; }
48
50 const G4double& userTimeStep = 0);
51
52 void StartTracking(G4Track* aTrack) override;
53
55 const G4Track&, // track
56 G4double, // previousStepSize
57 G4ForceCondition* pForceCond) override;
58 G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&) override;
59
61 {
62 return DBL_MAX;
63 }
64
67 inline void SetVerbose(G4int verbose);
68
69 protected:
71
72 private:
73 struct G4PolyNucleotideReactionState
74 : public G4ProcessStateBase<G4PolyNucleotideReactionState>
75 {
76 G4PolyNucleotideReactionState();
77 ~G4PolyNucleotideReactionState() override = default;
78 G4String GetType() override { return "PolyNucleotideReactionState"; }
79
80 using DNANode =
81 std::variant<const G4DNAComponentNode*, /*for dnadamage chain*/
82 const G4VPhysicalVolume* /*for molecularDNA chain*/>;
83 DNANode fNodeReactant;
84 G4double fSampledMinTimeStep;
85 G4double fPreviousTimeAtPreStepPoint;
86 };
87 G4bool fHasAlreadyReachedNullTime;
88 G4int fVerbose;
89 G4double fRCutOff;
90 G4VDNAHitModel* fpDamageModel;
92};
94{
95 fVerbose = verbose;
96}
97
99 G4VDNAHitModel* pModel)
100{
101 fpDamageModel = pModel;
102}
103#endif
G4ForceCondition
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4DNAPolyNucleotideReactionProcess & operator=(const G4DNAPolyNucleotideReactionProcess &)=delete
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *pForceCond) override
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *) override
G4double CalculateTimeStep(const G4Track &trackA, const G4double &userTimeStep=0)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &) override
G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: G4Step.hh:62
#define DBL_MAX
Definition: templates.hh:62