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
parameterisation
include
G4VFastSimulationModel.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
// $Id$
28
//
29
//
30
//---------------------------------------------------------------
31
//
32
// G4VFastSimulationModel.hh
33
//
34
// Description:
35
// Base class for fast simulation models.
36
//
37
// History:
38
// Oct 97: Verderi && MoraDeFreitas - First Implementation.
39
//
40
//---------------------------------------------------------------
41
42
43
#ifndef G4VFastSimulationModel_h
44
#define G4VFastSimulationModel_h
45
46
#include "
G4FastTrack.hh
"
47
#include "
G4FastStep.hh
"
48
49
//---------------------------
50
// For possible future needs:
51
//---------------------------
52
typedef
G4Region
G4Envelope
;
53
54
//-------------------------------------------
55
//
56
// G4VFastSimulationModel class
57
//
58
//-------------------------------------------
59
60
// Class Description:
61
// This is the abstract class for the implementation of parameterisations.
62
// You have to inherit from it to implement your concrete parameterisation
63
// model.
64
//
65
66
class
G4VFastSimulationModel
67
{
68
public
:
// With description
69
70
G4VFastSimulationModel
(
const
G4String
& aName);
71
// aName identifies the parameterisation model.
72
73
G4VFastSimulationModel
(
const
G4String
& aName,
G4Envelope
*,
74
G4bool
IsUnique=
FALSE
);
75
// This constructor allows you to get a quick "getting started".
76
// In addition to the model name, this constructor accepts a G4LogicalVolume
77
// pointer. This volume will automatically becomes the envelope, and the
78
// needed G4FastSimulationManager object is constructed if necessary giving
79
// it the G4LogicalVolume pointer and the boolean value. If it already
80
// exists, the model is simply added to this manager. However the
81
// G4VFastSimulationModel object will not keep track of the envelope given
82
// in the constructor.
83
// The boolean argument is there for optimization purpose: if you know that
84
// the G4LogicalVolume envelope is placed only once you can turn this
85
// boolean value to "true" (an automated mechanism is foreseen here.)
86
87
public
:
// Without description
88
virtual
~G4VFastSimulationModel
() {};
89
90
public
:
// With description
91
92
virtual
G4bool
IsApplicable
(
const
G4ParticleDefinition
&) = 0;
93
// In your implementation, you have to return "true" when your model is
94
// applicable to the G4ParticleDefinition passed to this method. The
95
// G4ParticleDefinition provides all intrisic particle informations (mass,
96
// charge, spin, name ...).
97
98
virtual
G4bool
ModelTrigger
(
const
G4FastTrack
&) = 0;
99
// You have to return "true" when the dynamics conditions to trigger your
100
// parameterisation are fulfiled. The G4FastTrack provides you access to
101
// the current G4Track, gives simple access to envelope related features
102
// (G4LogicalVolume, G4VSolid, G4AffineTransform references between the
103
// global and the envelope local coordinates systems) and simple access to
104
// the position, momentum expressed in the envelope coordinate system.
105
// Using those quantities and the G4VSolid methods, you can for example
106
// easily check how far you are from the envelope boundary.
107
108
virtual
void
DoIt
(
const
G4FastTrack
&,
G4FastStep
&) = 0;
109
// Your parameterisation properly said. The G4FastTrack reference provides
110
// input informations. The final state of the particles after parameterisation
111
// has to be returned through the G4FastStep reference. This final state is
112
// described has "requests" the tracking will apply after your
113
// parameterisation has been invoked.
114
115
// ---------------------------
116
// -- Idem for AtRest methods:
117
// ---------------------------
118
// -- A default dummy implementation is provided.
119
120
virtual
121
G4bool
AtRestModelTrigger
(
const
G4FastTrack
&) {
return
false
;}
122
// You have to return "true" when the dynamics conditions to trigger your
123
// parameterisation are fulfiled. The G4FastTrack provides you access to
124
// the current G4Track, gives simple access to envelope related features
125
// (G4LogicalVolume, G4VSolid, G4AffineTransform references between the
126
// global and the envelope local coordinates systems) and simple access to
127
// the position, momentum expressed in the envelope coordinate system.
128
// Using those quantities and the G4VSolid methods, you can for example
129
// easily check how far you are from the envelope boundary.
130
131
virtual
132
void
AtRestDoIt
(
const
G4FastTrack
&,
G4FastStep
&) {}
133
// Your parameterisation properly said. The G4FastTrack reference provides
134
// input informations. The final state of the particles after parameterisation
135
// has to be returned through the G4FastStep reference. This final state is
136
// described has "requests" the tracking will apply after your
137
// parameterisation has been invoked.
138
139
public
:
// Without description
140
141
// Useful public methods :
142
const
G4String
GetName
()
const
;
143
G4bool
operator ==
(
const
G4VFastSimulationModel
&)
const
;
144
145
private
:
146
//-------------
147
// Model Name:
148
//-------------
149
G4String
theModelName;
150
};
151
152
inline
const
G4String
G4VFastSimulationModel::GetName
()
const
153
{
154
return
theModelName;
155
}
156
157
inline
G4bool
158
G4VFastSimulationModel::operator ==
(
const
G4VFastSimulationModel
& fsm)
const
159
{
160
return
(
this
==&fsm) ?
true
:
false
;
161
}
162
#endif
Generated on Sat May 25 2013 14:34:11 for Geant4 by
1.8.4