Geant4-11
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// ------------------------------------------
66
67// --------------------------------------------------
68// -- static methods to retrieve the manager pointer:
69// --------------------------------------------------
71{
74
76}
77
78
80{
82}
83
84// ---------------
85// -- constructor
86// ---------------
88{
90}
91
92// -------------
93// -- destructor
94// -------------
96{
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{
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
150
151// ---------------------------------
152// -- display fast simulation setup:
153// ---------------------------------
155{
156 std::vector<G4VPhysicalVolume*> worldDone;
157 G4VPhysicalVolume* world;
159 // ----------------------------------------------------
160 // -- loop on regions to get the list of world volumes:
161 // ----------------------------------------------------
162 G4cout << "\nFast simulation setup:" << G4endl;
163 for (size_t i=0; i<regions->size(); i++)
164 {
165 world = (*regions)[i]->GetWorldPhysical();
166 if (world == nullptr) // region does not belong to any (existing) world
167 {
168 continue;
169 }
170 G4bool newWorld = true;
171 for (size_t ii=0; ii<worldDone.size(); ii++) if (worldDone[ii] == world) {newWorld = false; break;}
172 if (newWorld)
173 {
174 worldDone.push_back(world);
175 G4Region* worldRegion = world->GetLogicalVolume()->GetRegion();
176 // -- preambule: print physical volume and region names...
177 if (world == G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume())
178 G4cout << "\n * Mass Geometry with ";
179 else
180 G4cout << "\n * Parallel Geometry with ";
181 G4cout << "world volume: `" << world->GetName() << "' [region : `" << worldRegion->GetName() << "']" << G4endl;
182 // -- ... and print G4FSMP(s) attached to this world volume:
183 G4bool findG4FSMP(false);
184 // -- show to what particles this G4FSMP is attached to:
185 std::vector<G4ParticleDefinition*> particlesKnown;
186 for (size_t ip=0; ip<fFSMPVector.size(); ip++)
187 if (fFSMPVector[ip]->GetWorldVolume() == world)
188 {
189 G4cout << " o G4FastSimulationProcess: '" << fFSMPVector[ip]->GetProcessName() << "'" << G4endl;
190 G4cout << " Attached to:";
192 for (G4int iParticle=0; iParticle<particles->entries(); iParticle++)
193 {
194 G4ParticleDefinition* particle = particles->GetParticle(iParticle);
195 G4ProcessVector* processes = particle->GetProcessManager()->GetProcessList();
196 if (processes->contains(fFSMPVector[ip])) {G4cout << " " << particle->GetParticleName(); findG4FSMP = true; particlesKnown.push_back(particle);}
197 }
198 G4cout << G4endl;
199 }
200 if (!findG4FSMP) G4cout << " o G4FastSimulationProcess: (none)" << G4endl;
201 // -- now display the regions in this world volume, with mother<->daughter link shown by indentation:
202 G4cout << " o Region(s) and model(s) setup:" << G4endl;
203 DisplayRegion(worldRegion, 1, particlesKnown);
204 }
205 }
206}
207
208
209void G4GlobalFastSimulationManager::DisplayRegion(G4Region* region, G4int depth, std::vector<G4ParticleDefinition*>& particlesKnown) const
210{
211 G4String indent = " ";
212 for (G4int I=0; I<depth; I++) indent += " ";
213 G4cout << indent << "Region: `" << region->GetName() <<"'" << G4endl;
214 G4FastSimulationManager* fastSimManager = region->GetFastSimulationManager();
215 if (fastSimManager)
216 {
217 indent += " ";
218 G4cout << indent << "Model(s):" << G4endl;
219 indent += " ";
220 for (size_t im=0; im<fastSimManager->GetFastSimulationModelList().size(); im++)
221 {
222 G4cout << indent << "`" << (fastSimManager->GetFastSimulationModelList())[im]->GetName() << "'";
223 G4cout << " ; applicable to:";
225 for (G4int iParticle=0; iParticle<particles->entries(); iParticle++)
226 {
227 if ((fastSimManager->GetFastSimulationModelList())[im]->IsApplicable(*(particles->GetParticle(iParticle))))
228 {
229 G4cout << " " << particles->GetParticle(iParticle)->GetParticleName();
230 G4bool known(false);
231 for (size_t l=0; l<particlesKnown.size();l++) if(particlesKnown[l] == particles->GetParticle(iParticle)) {known = true; break;}
232 if (!known) G4cout << "[!!]";
233 }
234 }
235 G4cout << G4endl;
236 }
237 }
238
239 // -- all that to check mothership of "region"
241 for (size_t ip=0; ip<physVolStore->size(); ip++)
242 {
243 G4VPhysicalVolume* physVol = (*physVolStore)[ip];
244 if (physVol->GetLogicalVolume()->IsRootRegion())
245 if (physVol->GetMotherLogical())
246 {
247 G4Region* thisVolMotherRegion = physVol->GetMotherLogical()->GetRegion();
248 if (thisVolMotherRegion == region)
249 DisplayRegion(physVol->GetLogicalVolume()->GetRegion(), depth+1, particlesKnown);
250 }
251 }
252}
253
254
255// ----------------------------
256// -- management methods : list
257// ----------------------------
258
260 listType theType)
261{
262 if (theType == ISAPPLICABLE)
263 {
264 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++) ManagedManagers[ifsm]->ListModels(aName);
265 return;
266 }
267
268 if(aName == "all")
269 {
270 G4int titled = 0;
271 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
272 {
273 if(theType == NAMES_ONLY)
274 {
275 if(!(titled++))
276 G4cout << "Current Envelopes for Fast Simulation:\n";
277 G4cout << " ";
278 ManagedManagers[ifsm]->ListTitle();
279 G4cout << G4endl;
280 }
281 else ManagedManagers[ifsm]->ListModels();
282 }
283 }
284 else
285 {
286 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
287 if(aName == ManagedManagers[ifsm]-> GetEnvelope()->GetName())
288 {
289 ManagedManagers[ifsm]->ListModels();
290 break;
291 }
292 }
293}
294
296{
297 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
298 ManagedManagers[ifsm]->ListModels(aPD);
299}
300
301
303 const G4VFastSimulationModel* previousFound) const
304{
305 G4VFastSimulationModel* model = 0;
306 // -- flag used to navigate accross the various managers;
307 bool foundPrevious(false);
308 for (size_t ifsm=0; ifsm<ManagedManagers.size(); ifsm++)
309 {
310 model = ManagedManagers[ifsm]->
311 GetFastSimulationModel(modelName, previousFound, foundPrevious);
312 if (model) break;
313 }
314 return model;
315}
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 *)
G4FastSimulationVector< G4FastSimulationManagerProcess > fFSMPVector
void RemoveFSMP(G4FastSimulationManagerProcess *)
G4VFastSimulationModel * GetFastSimulationModel(const G4String &modelName, const G4VFastSimulationModel *previousFound=0) const
void ListEnvelopes(const G4String &aName="all", listType aListType=NAMES_ONLY)
void DisplayRegion(G4Region *motherRegion, G4int depth, std::vector< G4ParticleDefinition * > &particles) const
void InActivateFastSimulationModel(const G4String &)
G4FastSimulationVector< G4FastSimulationManager > ManagedManagers
void AddFSMP(G4FastSimulationManagerProcess *)
G4FastSimulationMessenger * fTheFastSimulationMessenger
void RemoveFastSimulationManager(G4FastSimulationManager *)
static G4ThreadLocal G4GlobalFastSimulationManager * fGlobalFastSimulationManager
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