source: trunk/source/processes/parameterisation/include/G4VFastSimulationModel.hh@ 1066

Last change on this file since 1066 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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