Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNADamage.cc
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#include "G4DNADamage.hh"
28#include "G4UnitsTable.hh"
29
31
33 const G4Molecule* molecule,
35 G4double time) :
36 G4VDNAHit(), fpMolecule(molecule)
37{
38 fBaseName = baseName;
40 fTime = time;
41}
42
44{
45 if (fpMolecule) delete fpMolecule;
46 fpMolecule = 0;
47}
48
50{
51 G4cout << "Reaction : " << fpMolecule->GetName() << " + " << fBaseName
52 << " at position : " << G4BestUnit(fPosition, "Length")
53 << " and time : " << G4BestUnit(fTime, "Time") << G4endl;
54}
55
57{
58 if (!fpInstance) new G4DNADamage();
59
60 return fpInstance;
61}
62
64{
65 fJustCountDamage = false;
67 fpInstance = this;
68}
69
71{
72 for (int i = 0; i < (int) fIndirectHits.size(); i++)
73 {
74 if (fIndirectHits[i]) delete fIndirectHits[i];
75 }
76 fIndirectHits.clear();
77}
78
80{
81 if (fpInstance) delete fpInstance;
82 fpInstance = 0;
83}
84
86{
88 for (int i = 0; i < (int) fIndirectHits.size(); i++)
89 {
90 if (fIndirectHits[i]) delete fIndirectHits[i];
91 }
92 fIndirectHits.clear();
93}
94
96 const G4Molecule* molecule,
98 G4double time)
99{
101 {
103 return;
104 }
105
106 G4DNAIndirectHit* indirectHit = 0;
107 std::map<G4Molecule, const G4Molecule*>::iterator it =
108 fMolMap.find(*molecule);
109
110 if (it == fMolMap.end())
111 {
112 G4Molecule* mol(0);
113 fMolMap[*molecule] = (mol = new G4Molecule(*molecule));
114 indirectHit = new G4DNAIndirectHit(baseName, mol, position, time);
115 }
116 else
117 {
118 indirectHit = new G4DNAIndirectHit(baseName, it->second, position, time);
119 }
120 fIndirectHits.push_back(indirectHit);
121}
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void Reset()
Definition: G4DNADamage.cc:85
virtual ~G4DNADamage()
Definition: G4DNADamage.cc:70
std::map< G4Molecule, const G4Molecule * > fMolMap
Definition: G4DNADamage.hh:120
G4int fNIndirectDamage
Definition: G4DNADamage.hh:118
virtual void AddIndirectDamage(const G4String &baseName, const G4Molecule *molecule, const G4ThreeVector &position, double time)
Definition: G4DNADamage.cc:95
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
G4DNAIndirectHit(const G4String &baseName, const G4Molecule *molecule, const G4ThreeVector &position, G4double time)
Definition: G4DNADamage.cc:32
virtual ~G4DNAIndirectHit()
Definition: G4DNADamage.cc:43
G4String fBaseName
Definition: G4DNADamage.hh:76
const G4Molecule * fpMolecule
Definition: G4DNADamage.hh:73
G4ThreeVector fPosition
Definition: G4DNADamage.hh:74
const G4String & GetName() const
Definition: G4Molecule.cc:338
#define G4ThreadLocal
Definition: tls.hh:77
#define position
Definition: xmlparse.cc:622