Geant4-11
G4FastSimulationManager.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//
28//
29//---------------------------------------------------------------
30//
31// G4FastSimulationManager.hh
32//
33// Description:
34// Manages the Fast Simulation models attached to a envelope.
35//
36// History:
37// Oct 97: Verderi && MoraDeFreitas - First Implementation.
38//
39//---------------------------------------------------------------
40
41
42#ifndef G4FastSimulationManager_h
43#define G4FastSimulationManager_h 1
44
45#include "globals.hh"
46
47#include "G4LogicalVolume.hh"
48#include "G4Region.hh"
49#include "G4VPhysicalVolume.hh"
50#include "G4ParticleTable.hh"
52#include "G4VParticleChange.hh"
53#include "G4FastTrack.hh"
54#include "G4FastStep.hh"
56#include "G4RotationMatrix.hh"
57#include "G4ThreeVector.hh"
58#include "G4Transform3D.hh"
60
61#include "G4ios.hh"
62
63//-------------------------------------------
64//
65// G4FastSimulationManager class
66//
67//-------------------------------------------
68
69// Class Description:
70// The G4VFastSimulationModel objects are attached to the envelope
71// through a G4FastSimulationManager.
72// This object will manage the list of models and will message them
73// at tracking time.
74//
75
76
78{
79public: // with description
80 //------------------------
81 // Constructor/Destructor
82 //------------------------
83 // Only one Constructor. By default the envelope can
84 // be placed n-Times.
85 // If the user is sure that it is placed just one time,
86 // the IsUnique flag should be set TRUE to avoid the
87 // G4AffineTransform re-calculations each time we reach
88 // the envelope.
89
91 G4bool IsUnique = FALSE);
92 // This is the only constructor. In this constructor you specify
93 // the envelope by giving the G4Region (typedef-ed as G4Envelope)
94 // pointer. The G4FastSimulationManager object will bind itself to
95 // this envelope and will notify this G4Region to become an envelope.
96 // If you know that this region is used for only one logical volume,
97 // you can turn the IsUnique boolean to "true" to allow some optimization.
98 //
99 // Note that if you choose to use the G4VFastSimulationModel(const G4String&,
100 // G4Region*, G4bool) constructor for you model, the G4FastSimulationManager
101 // will be constructed using the given G4Region* and G4bool values of the
102 // model constructor.
103 //
104
105public: // without description
107
108
109public: // with description
110 // Methods to add/remove models to/from the Model
111 // List.
112 //
114 // Add a model to the Model List.
115
117 // Remove a model from the Model List.
118
119 // Methods to activate/inactivate models from the Model
120 // List.
121
123 // Activate a model in the Model List.
124
126 // Inactivate a model in the Model List.
127
128public: // without description
129 // Methods for print/control commands
130 void ListTitle() const;
131 void ListModels() const;
132 void ListModels(const G4ParticleDefinition*) const;
133 void ListModels(const G4String& aName) const;
134 const G4Envelope* GetEnvelope() const;
135
137 const G4VFastSimulationModel* previousFound,
138 bool &foundPrevious) const;
139
140 const std::vector<G4VFastSimulationModel*>& GetFastSimulationModelList() const
141 {return ModelList;}
142
143
144 //----------------------------------------------
145 // Interface methods for the
146 // G4FastSimulationManagerProcess process.
147 //----------------------------------------------
148 // Trigger
150 const G4Navigator* a = 0);
151 // DoIt
153
154 // AtRest methods:
156 const G4Navigator* a = 0);
158
159 // For management
161
162private:
163 // Private members :
169
172
173 // -- *** depracating, to be dropped @ next major release:
175};
176
177inline void
179{
180 ModelList.push_back(fsm);
181 // forces the fApplicableModelList to be rebuild
183}
184
185inline void
187{
189 // forces the fApplicableModelList to be rebuild
191}
192
193inline G4bool
195{
196 return (this==&fsm) ? true : false;
197}
198
199inline const G4Envelope*
201{
202 return fFastTrack.GetEnvelope();
203}
204
205#endif
bool G4bool
Definition: G4Types.hh:86
#define FALSE
Definition: Globals.hh:23
G4FastSimulationManager(G4Envelope *anEnvelope, G4bool IsUnique=FALSE)
void RemoveFastSimulationModel(G4VFastSimulationModel *)
const G4Envelope * GetEnvelope() const
G4bool operator==(const G4FastSimulationManager &) const
G4VParticleChange * InvokePostStepDoIt()
G4FastSimulationVector< G4VFastSimulationModel > fInactivatedModels
G4FastSimulationVector< G4VFastSimulationModel > fApplicableModelList
G4VParticleChange * InvokeAtRestDoIt()
const std::vector< G4VFastSimulationModel * > & GetFastSimulationModelList() const
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4VFastSimulationModel * fTriggedFastSimulationModel
G4bool ActivateFastSimulationModel(const G4String &)
G4VFastSimulationModel * GetFastSimulationModel(const G4String &modelName, const G4VFastSimulationModel *previousFound, bool &foundPrevious) const
G4FastSimulationVector< G4VFastSimulationModel > ModelList
G4FastSimulationVector< G4Transform3D > GhostPlacements
G4bool InActivateFastSimulationModel(const G4String &)
void AddFastSimulationModel(G4VFastSimulationModel *)
G4ParticleDefinition * fLastCrossedParticle
T * remove(const T *)
G4Envelope * GetEnvelope() const
Definition: G4FastTrack.hh:186