Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
dna
processes
include
G4DNAMolecularDecay.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: G4DNAMolecularDecay.hh 65022 2012-11-12 16:43:12Z gcosmo $
27
//
28
// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29
//
30
// WARNING : This class is released as a prototype.
31
// It might strongly evolve or even disapear in the next releases.
32
//
33
// History:
34
// -----------
35
// 10 Oct 2011 M.Karamitros created
36
//
37
// -------------------------------------------------------------------
38
39
#ifndef G4MOLECULARDECAYPROCESS_HH
40
#define G4MOLECULARDECAYPROCESS_HH
41
42
#include "
G4VITRestProcess.hh
"
43
#include <map>
44
45
class
G4ParticleChange
;
46
class
G4VMolecularDecayDisplacer
;
47
class
G4MoleculeDefinition
;
48
56
class
G4DNAMolecularDecay
:
public
G4VITRestProcess
57
{
58
public
:
59
G4DNAMolecularDecay
(
const
G4String
& processName =
"DNAMolecularDecay"
,
60
G4ProcessType
type =
fDecay
);
61
62
virtual
~G4DNAMolecularDecay
();
63
64
G4IT_ADD_CLONE
(
G4VITProcess
,
G4DNAMolecularDecay
)
65
66
virtual
G4bool
IsApplicable
(
const
G4ParticleDefinition
&);
67
68
inline
G4double
AtRestGetPhysicalInteractionLength
(
69
const
G4Track
& track,
70
G4ForceCondition
*
condition
71
);
72
inline
G4VParticleChange
*
AtRestDoIt
(
73
const
G4Track
& track,
74
const
G4Step
& step
75
);
76
77
inline
void
SetVerbose
(
G4int
);
78
79
//__________________________________________________________________
80
void
SetDecayDisplacer
(
const
G4ParticleDefinition
*,
G4VMolecularDecayDisplacer
*);
81
G4VMolecularDecayDisplacer
*
GetDecayDisplacer
(
const
G4ParticleDefinition
*);
82
83
protected
:
84
//__________________________________________________________________
85
// Make the decay
86
virtual
G4VParticleChange
*
DecayIt
(
const
G4Track
& ,
const
G4Step
&);
87
virtual
G4double
GetMeanLifeTime
(
const
G4Track
&,
G4ForceCondition
*);
88
89
private
:
90
G4DNAMolecularDecay
();
91
G4DNAMolecularDecay
(
const
G4DNAMolecularDecay
&
right
);
92
G4DNAMolecularDecay
& operator=(
const
G4DNAMolecularDecay
&right);
93
94
private
:
95
G4bool
fDecayAtFixedTime ;
96
97
typedef
std::map<const G4ParticleDefinition*, G4VMolecularDecayDisplacer*> DecayDisplacementMap;
98
DecayDisplacementMap fDecayDisplacementMap;
99
100
G4int
fVerbose;
101
};
102
103
inline
void
G4DNAMolecularDecay::SetVerbose
(
G4int
verbose)
104
{
105
fVerbose = verbose ;
106
}
107
108
inline
G4double
G4DNAMolecularDecay::AtRestGetPhysicalInteractionLength
(
109
const
G4Track
& track,
110
G4ForceCondition
*
condition
)
111
{
112
if
(fDecayAtFixedTime)
113
{
114
return
GetMeanLifeTime
(track, condition);
115
}
116
117
return
G4VITRestProcess::AtRestGetPhysicalInteractionLength
(track, condition);
118
}
119
120
inline
G4VParticleChange
*
G4DNAMolecularDecay::AtRestDoIt
(
121
const
G4Track
& track,
122
const
G4Step
& step
123
)
124
{
125
ClearNumberOfInteractionLengthLeft
();
126
ClearInteractionTimeLeft
();
127
return
DecayIt
(track, step);
128
}
129
130
131
#endif
/* G4MOLECULARDECAYPROCESS_HH */
Generated on Sat May 25 2013 14:33:31 for Geant4 by
1.8.4