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