Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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#ifndef G4FastSimulationManager_h
42#define G4FastSimulationManager_h 1
43
45#include "G4FastStep.hh"
46#include "G4FastTrack.hh"
47#include "G4LogicalVolume.hh"
49#include "G4ParticleTable.hh"
50#include "G4Region.hh"
51#include "G4RotationMatrix.hh"
52#include "G4ThreeVector.hh"
53#include "G4Transform3D.hh"
55#include "G4VParticleChange.hh"
56#include "G4VPhysicalVolume.hh"
57#include "G4ios.hh"
58#include "globals.hh"
59
60//-------------------------------------------
61//
62// G4FastSimulationManager class
63//
64//-------------------------------------------
65
66// Class Description:
67// The G4VFastSimulationModel objects are attached to the envelope
68// through a G4FastSimulationManager.
69// This object will manage the list of models and will message them
70// at tracking time.
71//
72
74{
75 public: // with description
76 //------------------------
77 // Constructor/Destructor
78 //------------------------
79 // Only one Constructor. By default the envelope can
80 // be placed n-Times.
81 // If the user is sure that it is placed just one time,
82 // the IsUnique flag should be set TRUE to avoid the
83 // G4AffineTransform re-calculations each time we reach
84 // the envelope.
85
86 G4FastSimulationManager(G4Envelope* anEnvelope, G4bool IsUnique = FALSE);
87 // This is the only constructor. In this constructor you specify
88 // the envelope by giving the G4Region (typedef-ed as G4Envelope)
89 // pointer. The G4FastSimulationManager object will bind itself to
90 // this envelope and will notify this G4Region to become an envelope.
91 // If you know that this region is used for only one logical volume,
92 // you can turn the IsUnique boolean to "true" to allow some optimization.
93 //
94 // Note that if you choose to use the G4VFastSimulationModel(const G4String&,
95 // G4Region*, G4bool) constructor for you model, the G4FastSimulationManager
96 // will be constructed using the given G4Region* and G4bool values of the
97 // model constructor.
98 //
99
100 // Destructor
102
103 // Add a model to the Model List.
105
106 // Remove a model from the Model List.
108
109 // Activate a model in the Model List.
111
112 // Inactivate a model in the Model List.
114
115 // Methods for print/control commands
116 void ListTitle() const;
117 void ListModels() const;
118 void ListModels(const G4ParticleDefinition*) const;
119 void ListModels(const G4String& aName) const;
120 const G4Envelope* GetEnvelope() const;
121
123 const G4VFastSimulationModel* previousFound,
124 G4bool& foundPrevious) const;
125
126 const std::vector<G4VFastSimulationModel*>& GetFastSimulationModelList() const
127 {
128 return ModelList;
129 }
130
131 void FlushModels();
132
133 //----------------------------------------------
134 // Interface methods for the
135 // G4FastSimulationManagerProcess process.
136 //----------------------------------------------
137 // Trigger
139 // DoIt
141
142 // AtRest methods:
145
146 // For management
148
149 private:
150 // Private members :
151 G4FastTrack fFastTrack;
152 G4FastStep fFastStep;
153 G4VFastSimulationModel* fTriggedFastSimulationModel{nullptr};
156
157 G4ParticleDefinition* fLastCrossedParticle{nullptr};
159
160 // -- *** depracating, to be dropped @ next major release:
162};
163
165{
166 ModelList.push_back(fsm);
167 // forces the fApplicableModelList to be rebuild
168 fLastCrossedParticle = nullptr;
169}
170
172{
173 if (ModelList.remove(fsm) == nullptr) fInactivatedModels.remove(fsm);
174 // forces the fApplicableModelList to be rebuild
175 fLastCrossedParticle = nullptr;
176}
177
179{
180 return this == &fsm;
181}
182
184{
185 return fFastTrack.GetEnvelope();
186}
187
188#endif
bool G4bool
Definition G4Types.hh:86
G4FastSimulationManager(G4Envelope *anEnvelope, G4bool IsUnique=FALSE)
void RemoveFastSimulationModel(G4VFastSimulationModel *)
const G4Envelope * GetEnvelope() const
G4bool operator==(const G4FastSimulationManager &) const
G4VParticleChange * InvokePostStepDoIt()
G4VParticleChange * InvokeAtRestDoIt()
const std::vector< G4VFastSimulationModel * > & GetFastSimulationModelList() const
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)
G4bool ActivateFastSimulationModel(const G4String &)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)
G4bool InActivateFastSimulationModel(const G4String &)
void AddFastSimulationModel(G4VFastSimulationModel *)
G4VFastSimulationModel * GetFastSimulationModel(const G4String &modelName, const G4VFastSimulationModel *previousFound, G4bool &foundPrevious) const
T * remove(const T *)
G4Envelope * GetEnvelope() const
#define FALSE
Definition globals.hh:38