Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GlobalFastSimulationManager.cc
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// G4GlobalFastSimulationManager.cc
32//
33// Description:
34// A singleton class which manages the Fast Simulation managers
35// attached to envelopes. Implementation.
36//
37// History:
38// June 98: Verderi && MoraDeFreitas - "G4ParallelWorld" becomes
39// "G4FlavoredParallelWorld"; some method name changes;
40// GetFlavoredWorldForThis now returns a
41// G4FlavoredParallelWorld pointer.
42// Feb 98: Verderi && MoraDeFreitas - First Implementation.
43// March 98: correction to instanciate dynamically the manager
44// May 07: Move to parallel world scheme
45//
46//---------------------------------------------------------------
47
49#include "G4ParticleTable.hh"
51#include "G4Material.hh"
52#include "G4ThreeVector.hh"
53#include "G4PVPlacement.hh"
56#include "G4RegionStore.hh"
57#include "G4ProcessVector.hh"
58#include "G4ProcessManager.hh"
60
61
62// ------------------------------------------
63// -- static instance pointer initialisation:
64// ------------------------------------------
65G4ThreadLocal G4GlobalFastSimulationManager* G4GlobalFastSimulationManager::fGlobalFastSimulationManager = 0;
66
67// --------------------------------------------------
68// -- static methods to retrieve the manager pointer:
69// --------------------------------------------------
71{
72 if(!fGlobalFastSimulationManager)
73 fGlobalFastSimulationManager = new G4GlobalFastSimulationManager;
74
75 return fGlobalFastSimulationManager;
76}
77
78
80{
82}
83
84// ---------------
85// -- constructor
86// ---------------
87G4GlobalFastSimulationManager::G4GlobalFastSimulationManager()
88{
89 fTheFastSimulationMessenger = new G4FastSimulationMessenger(this);
90}
91
92// -------------
93// -- destructor
94// -------------
96{
97 delete fTheFastSimulationMessenger;
98 fTheFastSimulationMessenger = 0;
99}
100
101// ----------------------
102// -- management methods:
103// ----------------------
106{
107 ManagedManagers.push_back(fsmanager);
108}
109
112{
113 ManagedManagers.remove(fsmanager);
114}
115
117{
118 fFSMPVector.push_back(fp);
119}
120
122{
123 fFSMPVector.remove(fp);
124}
125
127{
128 G4bool result = false;
129 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
130 result = result || ManagedManagers[ifsm]->
132 if(result)
133 G4cout << "Model " << aName << " activated.";
134 else
135 G4cout << "Model " << aName << " not found.";
136 G4cout << G4endl;
137}
138
140{
141 G4bool result = false;
142 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
143 result = result || ManagedManagers[ifsm]->
145 if (result) G4cout << "Model " << aName << " inactivated.";
146 else G4cout << "Model " << aName << " not found.";
147 G4cout << G4endl;
148}
149
151{
152 // loop over all models (that need flushing?) and flush
153 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
154 ManagedManagers[ifsm]->FlushModels();
155}
156
157// ---------------------------------
158// -- display fast simulation setup:
159// ---------------------------------
161{
162 std::vector<G4VPhysicalVolume*> worldDone;
163 G4VPhysicalVolume* world;
165 // ----------------------------------------------------
166 // -- loop on regions to get the list of world volumes:
167 // ----------------------------------------------------
168 G4cout << "\nFast simulation setup:" << G4endl;
169 for (size_t i=0; i<regions->size(); i++)
170 {
171 world = (*regions)[i]->GetWorldPhysical();
172 if (world == nullptr) // region does not belong to any (existing) world
173 {
174 continue;
175 }
176 G4bool newWorld = true;
177 for (size_t ii=0; ii<worldDone.size(); ii++) if (worldDone[ii] == world) {newWorld = false; break;}
178 if (newWorld)
179 {
180 worldDone.push_back(world);
181 G4Region* worldRegion = world->GetLogicalVolume()->GetRegion();
182 // -- preambule: print physical volume and region names...
183 if (world == G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume())
184 G4cout << "\n * Mass Geometry with ";
185 else
186 G4cout << "\n * Parallel Geometry with ";
187 G4cout << "world volume: `" << world->GetName() << "' [region : `" << worldRegion->GetName() << "']" << G4endl;
188 // -- ... and print G4FSMP(s) attached to this world volume:
189 G4bool findG4FSMP(false);
190 // -- show to what particles this G4FSMP is attached to:
191 std::vector<G4ParticleDefinition*> particlesKnown;
192 for (size_t ip=0; ip<fFSMPVector.size(); ip++)
193 if (fFSMPVector[ip]->GetWorldVolume() == world)
194 {
195 G4cout << " o G4FastSimulationProcess: '" << fFSMPVector[ip]->GetProcessName() << "'" << G4endl;
196 G4cout << " Attached to:";
198 for (G4int iParticle=0; iParticle<particles->entries(); iParticle++)
199 {
200 G4ParticleDefinition* particle = particles->GetParticle(iParticle);
201 G4ProcessVector* processes = particle->GetProcessManager()->GetProcessList();
202 if (processes->contains(fFSMPVector[ip])) {G4cout << " " << particle->GetParticleName(); findG4FSMP = true; particlesKnown.push_back(particle);}
203 }
204 G4cout << G4endl;
205 }
206 if (!findG4FSMP) G4cout << " o G4FastSimulationProcess: (none)" << G4endl;
207 // -- now display the regions in this world volume, with mother<->daughter link shown by indentation:
208 G4cout << " o Region(s) and model(s) setup:" << G4endl;
209 DisplayRegion(worldRegion, 1, particlesKnown);
210 }
211 }
212}
213
214
215void G4GlobalFastSimulationManager::DisplayRegion(G4Region* region, G4int depth, std::vector<G4ParticleDefinition*>& particlesKnown) const
216{
217 G4String indent = " ";
218 for (G4int I=0; I<depth; I++) indent += " ";
219 G4cout << indent << "Region: `" << region->GetName() <<"'" << G4endl;
220 G4FastSimulationManager* fastSimManager = region->GetFastSimulationManager();
221 if (fastSimManager)
222 {
223 indent += " ";
224 G4cout << indent << "Model(s):" << G4endl;
225 indent += " ";
226 for (size_t im=0; im<fastSimManager->GetFastSimulationModelList().size(); im++)
227 {
228 G4cout << indent << "`" << (fastSimManager->GetFastSimulationModelList())[im]->GetName() << "'";
229 G4cout << " ; applicable to:";
231 for (G4int iParticle=0; iParticle<particles->entries(); iParticle++)
232 {
233 if ((fastSimManager->GetFastSimulationModelList())[im]->IsApplicable(*(particles->GetParticle(iParticle))))
234 {
235 G4cout << " " << particles->GetParticle(iParticle)->GetParticleName();
236 G4bool known(false);
237 for (size_t l=0; l<particlesKnown.size();l++) if(particlesKnown[l] == particles->GetParticle(iParticle)) {known = true; break;}
238 if (!known) G4cout << "[!!]";
239 }
240 }
241 G4cout << G4endl;
242 }
243 }
244
245 // -- all that to check mothership of "region"
247 for (size_t ip=0; ip<physVolStore->size(); ip++)
248 {
249 G4VPhysicalVolume* physVol = (*physVolStore)[ip];
250 if (physVol->GetLogicalVolume()->IsRootRegion())
251 if (physVol->GetMotherLogical())
252 {
253 G4Region* thisVolMotherRegion = physVol->GetMotherLogical()->GetRegion();
254 if (thisVolMotherRegion == region)
255 DisplayRegion(physVol->GetLogicalVolume()->GetRegion(), depth+1, particlesKnown);
256 }
257 }
258}
259
260
261// ----------------------------
262// -- management methods : list
263// ----------------------------
264
266 listType theType)
267{
268 if (theType == ISAPPLICABLE)
269 {
270 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++) ManagedManagers[ifsm]->ListModels(aName);
271 return;
272 }
273
274 if(aName == "all")
275 {
276 G4int titled = 0;
277 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
278 {
279 if(theType == NAMES_ONLY)
280 {
281 if(!(titled++))
282 G4cout << "Current Envelopes for Fast Simulation:\n";
283 G4cout << " ";
284 ManagedManagers[ifsm]->ListTitle();
285 G4cout << G4endl;
286 }
287 else ManagedManagers[ifsm]->ListModels();
288 }
289 }
290 else
291 {
292 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
293 if(aName == ManagedManagers[ifsm]-> GetEnvelope()->GetName())
294 {
295 ManagedManagers[ifsm]->ListModels();
296 break;
297 }
298 }
299}
300
302{
303 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
304 ManagedManagers[ifsm]->ListModels(aPD);
305}
306
307
309 const G4VFastSimulationModel* previousFound) const
310{
311 G4VFastSimulationModel* model = 0;
312 // -- flag used to navigate accross the various managers;
313 bool foundPrevious(false);
314 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
315 {
316 model = ManagedManagers[ifsm]->
317 GetFastSimulationModel(modelName, previousFound, foundPrevious);
318 if (model) break;
319 }
320 return model;
321}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const std::vector< G4VFastSimulationModel * > & GetFastSimulationModelList() const
T * remove(const T *)
void RemoveFSMP(G4FastSimulationManagerProcess *)
G4VFastSimulationModel * GetFastSimulationModel(const G4String &modelName, const G4VFastSimulationModel *previousFound=0) const
void ListEnvelopes(const G4String &aName="all", listType aListType=NAMES_ONLY)
void InActivateFastSimulationModel(const G4String &)
void AddFSMP(G4FastSimulationManagerProcess *)
void RemoveFastSimulationManager(G4FastSimulationManager *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
static G4GlobalFastSimulationManager * GetInstance()
void AddFastSimulationManager(G4FastSimulationManager *)
G4bool IsRootRegion() const
G4Region * GetRegion() const
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4ParticleDefinition * GetParticle(G4int index) const
G4int entries() const
static G4ParticleTable * GetParticleTable()
static G4PhysicalVolumeStore * GetInstance()
G4ProcessVector * GetProcessList() const
G4bool contains(G4VProcess *aProcess) const
static G4RegionStore * GetInstance()
G4FastSimulationManager * GetFastSimulationManager() const
Definition: G4Region.cc:140
const G4String & GetName() const
static G4TransportationManager * GetTransportationManager()
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
#define G4ThreadLocal
Definition: tls.hh:77