Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNADamage.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// Author: Mathieu Karamitros
28
29// The code is developed in the framework of the ESA AO7146
30//
31// We would be very happy hearing from you, send us your feedback! :)
32//
33// In order for Geant4-DNA to be maintained and still open-source,
34// article citations are crucial.
35// If you use Geant4-DNA chemistry and you publish papers about your software,
36// in addition to the general paper on Geant4-DNA:
37//
38// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39//
40// we would be very happy if you could please also cite the following
41// reference papers on chemistry:
42//
43// J. Comput. Phys. 274 (2014) 841-882
44// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45
46#ifndef G4DNADAMAGE_HH
47#define G4DNADAMAGE_HH 1
48
49#include "G4Molecule.hh"
50
52{
53public :
55 virtual ~G4VDNAHit(){;}
56};
57
59{
60public :
61 G4DNAIndirectHit(const G4String& baseName, const G4Molecule* molecule,
62 const G4ThreeVector& position, G4double time);
63 virtual ~G4DNAIndirectHit();
64
65 inline const G4Molecule* GetMolecule() {return fpMolecule;}
66 inline const G4ThreeVector& GetPosition() {return fPosition;}
67 inline const G4String& GetBaseName() {return fBaseName;}
68 inline double GetTime() {return fTime;}
69
70 void Print();
71
72protected :
77};
78
80{
81public:
82 static G4DNADamage* Instance();
83 static void DeleteInstance();
84
85 virtual void Reset();
86
87 //void AddDirectDamage();
88 virtual void AddIndirectDamage(const G4String& baseName,
89 const G4Molecule* molecule,
90 const G4ThreeVector& position,
91 G4double time);
92
93 inline const std::vector<G4DNAIndirectHit*>* GetIndirectHits();
94 inline virtual G4int GetNIndirectHits() const
95 {
97 return fNIndirectDamage;
98
99 return (G4int)fIndirectHits.size();
100 }
101
102 inline virtual void SetOnlyCountDamage(G4bool flag = true)
103 {
104 fJustCountDamage = flag;
105 }
106
107 inline virtual G4bool OnlyCountDamage() const
108 {
109 return fJustCountDamage;
110 }
111
112protected :
113 G4DNADamage();
115 virtual ~G4DNADamage();
116
119 std::vector<G4DNAIndirectHit*> fIndirectHits;
120 std::map<G4Molecule, const G4Molecule*> fMolMap;
121};
122
123inline const std::vector<G4DNAIndirectHit*>* G4DNADamage::GetIndirectHits()
124{
125 return &fIndirectHits;
126}
127
128#endif // G4DNADAMAGES_HH
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void Reset()
Definition: G4DNADamage.cc:85
virtual ~G4DNADamage()
Definition: G4DNADamage.cc:70
virtual G4bool OnlyCountDamage() const
Definition: G4DNADamage.hh:107
const std::vector< G4DNAIndirectHit * > * GetIndirectHits()
Definition: G4DNADamage.hh:123
virtual void AddIndirectDamage(const G4String &baseName, const G4Molecule *molecule, const G4ThreeVector &position, G4double time)
Definition: G4DNADamage.cc:95
std::map< G4Molecule, const G4Molecule * > fMolMap
Definition: G4DNADamage.hh:120
G4int fNIndirectDamage
Definition: G4DNADamage.hh:118
virtual void SetOnlyCountDamage(G4bool flag=true)
Definition: G4DNADamage.hh:102
virtual G4int GetNIndirectHits() const
Definition: G4DNADamage.hh:94
std::vector< G4DNAIndirectHit * > fIndirectHits
Definition: G4DNADamage.hh:119
static G4DNADamage * Instance()
Definition: G4DNADamage.cc:56
static void DeleteInstance()
Definition: G4DNADamage.cc:79
G4bool fJustCountDamage
Definition: G4DNADamage.hh:117
static G4ThreadLocal G4DNADamage * fpInstance
Definition: G4DNADamage.hh:114
double GetTime()
Definition: G4DNADamage.hh:68
const G4ThreeVector & GetPosition()
Definition: G4DNADamage.hh:66
virtual ~G4DNAIndirectHit()
Definition: G4DNADamage.cc:43
G4String fBaseName
Definition: G4DNADamage.hh:76
const G4Molecule * GetMolecule()
Definition: G4DNADamage.hh:65
const G4String & GetBaseName()
Definition: G4DNADamage.hh:67
const G4Molecule * fpMolecule
Definition: G4DNADamage.hh:73
G4ThreeVector fPosition
Definition: G4DNADamage.hh:74
virtual ~G4VDNAHit()
Definition: G4DNADamage.hh:55
#define G4ThreadLocal
Definition: tls.hh:77