Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FastSimulationManager.cc
이 파일의 문서화 페이지로 가기
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: G4FastSimulationManager.cc 100945 2016-11-03 11:26:19Z gcosmo $
28 //
29 //---------------------------------------------------------------
30 //
31 // G4FastSimulationManager.cc
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 // May 07: Move to parallel world scheme
40 //
41 //---------------------------------------------------------------
42 
45 #include "G4PVPlacement.hh"
47 
48 // --------------------------------------------------
49 // Constructor with envelope and IsUnique flag :
50 // --------------------------------------------------
51 //
54  G4bool IsUnique) :
55  fFastTrack(anEnvelope,IsUnique),fTriggedFastSimulationModel(0),
56  fLastCrossedParticle(0)
57 {
58  // Communicates to the region that it becomes a
59  // envelope and with this fast simulation manager.
60  anEnvelope->SetFastSimulationManager(this);
61 
62  // Add itself to the GlobalFastSimulationManager
64  AddFastSimulationManager(this);
65 }
66 
67 // -----------
68 // Destructor:
69 // -----------
71 {
72  //
73  // Check out the Envelope about this pointer. If in use,
74  // resets the Logical Volume IsEnvelope flag to avoid clash.
75  //
78  // Remove itself from the GlobalFastSimulationManager
80  RemoveFastSimulationManager(this);
81 }
82 
83 // ---------------------------------------
84 // Methods to activate/inactivate models
85 //----------------------------------------
86 
87 G4bool
89 {
90  size_t iModel;
91 
92  // If the model is already active, do nothing.
93  for (iModel=0; iModel<ModelList.size(); iModel++)
94  if(ModelList[iModel]->GetName() == aName)
95  return true;
96 
97  // Look for in the fInactivatedModels list, if found push_back it back to
98  // the ModelList
99  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
100  if(fInactivatedModels[iModel]->GetName() == aName) {
101  ModelList.
102  push_back (fInactivatedModels.removeAt(iModel));
103  // forces the fApplicableModelList to be rebuild
105  return true;
106  }
107  return false;
108 }
109 
110 G4bool
112 {
113  // Look for in the ModelList, if found remove from it and keep the pointer
114  // on the fInactivatedModels list.
115  for (size_t iModel=0; iModel<ModelList.size(); iModel++)
116  if(ModelList[iModel]->GetName() == aName) {
118  push_back (ModelList.removeAt(iModel));
119  // forces the fApplicableModelList to be rebuild
121  return true;
122  }
123  return false;
124 }
125 
128  const G4VFastSimulationModel* previousFound,
129  bool &foundPrevious) const
130 {
132  for (size_t iModel=0; iModel<ModelList.size(); iModel++)
133  {
134  if(ModelList[iModel]->GetName() == modelName)
135  {
136  if (previousFound == 0)
137  {
138  model = ModelList[iModel];
139  break;
140  }
141  else
142  {
143  if (ModelList[iModel] == previousFound)
144  {
145  foundPrevious = true;
146  continue;
147  }
148  if (foundPrevious)
149  {
150  model = ModelList[iModel];
151  break;
152  }
153  }
154  }
155  }
156  return model;
157 }
158 
159 
160 //------------------------------------------------------------------
161 // Interface trigger method for the G4ParameterisationManagerProcess
162 //------------------------------------------------------------------
163 // G4bool GetFastSimulationManagerTrigger(const G4Track &);
164 //
165 // This method is used to interface the G4FastSimulationManagerProcess
166 // with the user Fast Simulation Models. It's called when the particle
167 // is inside the envelope.
168 //
169 // It :
170 //
171 // 1) initialises the private members (fFastTrack and so
172 // on);
173 // 2) loops on the IsApplicable() methods to find out the
174 // ones should be applied.
175 // 2) for these, loops on the ModelTrigger() methods to find out
176 // perhaps one that must be applied just now.
177 //
178 // If the a Fast Simulation Model is triggered then it returns
179 // true, false otherwise.
180 //
181 //-----------------------------------------------------------
182 G4bool
185  const G4Navigator* theNavigator)
186 {
187  size_t iModel;
188 
189  // If particle type changed re-build the fApplicableModelList.
190  if(fLastCrossedParticle!=track.GetDefinition()) {
192  fApplicableModelList.clear();
193  // If Model List is empty, do nothing !
194  if(ModelList.size()==0) return false;
195  for (iModel=0; iModel<ModelList.size(); iModel++)
196  if(ModelList[iModel]->IsApplicable(*(track.GetDefinition())))
197  fApplicableModelList.push_back (ModelList[iModel]);
198  }
199 
200  // If Applicable Model List is empty, do nothing !
201  if(fApplicableModelList.size()==0) return false;
202 
203  // -- Register current track
204  fFastTrack.SetCurrentTrack(track,theNavigator);
205 
206  // tests if particle are on the boundary and leaving,
207  // in this case do nothing !
208  if(fFastTrack.OnTheBoundaryButExiting()) return false;
209 
210  // Loops on the ModelTrigger() methods
211  for (iModel=0; iModel<fApplicableModelList.size(); iModel++)
212 
213  //---------------------------------------------------
214  // Asks the ModelTrigger method if it must be trigged now.
215  //---------------------------------------------------
216 
217  if(fApplicableModelList[iModel]->ModelTrigger(fFastTrack)) {
218  //--------------------------------------------------
219  // The model will be applied. Initializes the G4FastStep
220  // with the current state of the G4Track and
221  // same usefull parameters.
222  // In particular it does SetLocalEnergyDeposit(0.0).
223  //--------------------------------------------------
225 
226  // Keeps the FastSimulationModel pointer to call the
227  // DoIt() method.
229  return true;
230  }
231 
232  //--------------------------------------------
233  // Nobody asks to gain control, returns false
234  //--------------------------------------------
235  return false;
236 }
237 
239 {
240  // const G4FastTrack& parFastTrack=fFastTrack;
242  return &fFastStep;
243 }
244 
245 // -------------------------------------------------------------
246 // -- Mostly the same as above, in the case of AtRest particles:
247 // -------------------------------------------------------------
248 G4bool
250  const G4Navigator* theNavigator)
251 {
252  size_t iModel;
253 
254  // If particle type changed re-build the fApplicableModelList.
255  if(fLastCrossedParticle!=track.GetDefinition()) {
257  fApplicableModelList.clear();
258  // If Model List is empty, do nothing !
259  if(ModelList.size()==0) return false;
260  for (iModel=0; iModel<ModelList.size(); iModel++)
261  if(ModelList[iModel]->IsApplicable(*(track.GetDefinition())))
262  fApplicableModelList.push_back (ModelList[iModel]);
263  }
264 
265  // If Applicable Model List is empty, do nothing !
266  if(fApplicableModelList.size()==0) return false;
267 
268  // -- Register current track
269  fFastTrack.SetCurrentTrack(track,theNavigator);
270 
271  // -- (note: compared to the PostStepGetFastSimulationManagerTrigger,
272  // -- the test to see if the particle is on the boundary but leaving
273  // -- is irrelevant here)
274 
275  // Loops on the models to see if one of them wants to trigger:
276  for (iModel=0; iModel < fApplicableModelList.size(); iModel++)
277  if(fApplicableModelList[iModel]->AtRestModelTrigger(fFastTrack))
278  {
281  return true;
282  }
283 
284  //--------------------------------------------
285  // Nobody asks to gain control, returns false
286  //--------------------------------------------
287  return false;
288 }
289 
291 {
293  return &fFastStep;
294 }
295 
296 void
298 {
300  // if(GhostPlacements.size()!=0) G4cout << " (ghost)";
302  else G4cout << " (// geom.)";
303 
304 }
305 
306 void
308 {
309  size_t iModel;
310 
311  G4cout << "Current Models for the ";
312  ListTitle();
313  G4cout << " envelope:\n";
314 
315  for (iModel=0; iModel<ModelList.size(); iModel++)
316  G4cout << " " << ModelList[iModel]->GetName() << "\n";
317 
318  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
319  G4cout << " " << fInactivatedModels[iModel]->GetName()
320  << "(inactivated)\n";
321 }
322 
323 void G4FastSimulationManager::ListModels(const G4String& modelName) const
324 {
325  size_t iModel;
326  G4int titled = 0;
328 
329  // Active Models
330  for ( iModel=0; iModel<ModelList.size(); iModel++ )
331  if( ModelList[iModel]->GetName() == modelName || modelName == "all" )
332  {
333  if( !(titled++) )
334  {
335  G4cout << "In the envelope ";
336  ListTitle();
337  G4cout << ",\n";
338  }
339  G4cout << " the model " << ModelList[iModel]->GetName()
340  << " is applicable for :\n ";
341 
342  G4int list_started=0;
343  for ( G4int iParticle = 0; iParticle<theParticleTable->entries(); iParticle++)
344  if( ModelList[iModel] -> IsApplicable( *(theParticleTable->GetParticle(iParticle))) )
345  {
346  if(list_started++) G4cout << ", ";
347  G4cout << theParticleTable->
348  GetParticle(iParticle)->GetParticleName();
349  }
350  G4cout <<G4endl;
351  }
352 
353  // Inactive Models
354  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
355  if(fInactivatedModels[iModel]->GetName() == modelName || modelName == "all" )
356  {
357  if( !(titled++) )
358  {
359  G4cout << "In the envelope ";
360  ListTitle();
361  G4cout << ",\n";
362  }
363  G4cout << " the model " << fInactivatedModels[iModel]->GetName()
364  << " (inactivated) is applicable for :\n ";
365 
366  G4int list_started=0;
367  for ( G4int iParticle=0; iParticle<theParticleTable->entries(); iParticle++ )
368  if( fInactivatedModels[iModel] -> IsApplicable( *(theParticleTable->GetParticle(iParticle))) )
369  {
370  if(list_started++) G4cout << ", ";
371  G4cout << theParticleTable->
372  GetParticle(iParticle)->GetParticleName();
373  }
374  G4cout <<G4endl;
375  }
376 }
377 
378 void G4FastSimulationManager::ListModels(const G4ParticleDefinition* particleDefinition) const
379 {
380  size_t iModel;
381  G4bool unique = true;
382 
383  // Active Models
384  for ( iModel=0; iModel<ModelList.size(); iModel++ )
385  if ( ModelList[iModel]->IsApplicable(*particleDefinition) )
386  {
387  G4cout << "Envelope ";
388  ListTitle();
389  G4cout << ", Model "
390  << ModelList[iModel]->GetName()
391  << "." << G4endl;
392  // -- Verify unicity of model attached to particleDefinition:
393  for ( auto jModel = iModel + 1; jModel < ModelList.size(); jModel++ )
394  if ( ModelList[jModel]->IsApplicable(*particleDefinition) ) unique = false;
395  }
396 
397  // Inactive Models
398  for ( iModel=0; iModel<fInactivatedModels.size(); iModel++ )
399  if( fInactivatedModels[iModel]->IsApplicable(*particleDefinition) )
400  {
401  G4cout << "Envelope ";
402  ListTitle();
403  G4cout << ", Model "
404  << fInactivatedModels[iModel]->GetName()
405  << " (inactivated)." << G4endl;
406  }
407 
408  if( !unique )
409  {
411  ed << "Two or more active Models are available for the same particle type, in the same envelope/region." << G4endl;
412  G4Exception("G4FastSimulationManager::ListModels(const G4ParticleDefinition* particleDefinition) const",
413  "FastSim001",
414  JustWarning, ed,
415  "Models risk to exclude each other.");
416  }
417  unique=false;
418 }
G4VParticleChange * InvokeAtRestDoIt()
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
virtual void DoIt(const G4FastTrack &, G4FastStep &)=0
G4FastSimulationManager * GetFastSimulationManager() const
Definition: G4Region.cc:137
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
static G4ParticleTable * GetParticleTable()
G4bool OnTheBoundaryButExiting() const
Definition: G4FastTrack.hh:243
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61
G4VParticleChange * InvokePostStepDoIt()
G4Envelope * GetEnvelope() const
Definition: G4FastTrack.hh:188
G4Navigator * GetNavigatorForTracking() const
G4FastSimulationVector< G4VFastSimulationModel > ModelList
G4FastSimulationVector< G4VFastSimulationModel > fApplicableModelList
G4ParticleDefinition * GetDefinition() const
G4bool ActivateFastSimulationModel(const G4String &)
G4int entries() const
void Initialize(const G4FastTrack &)
Definition: G4FastStep.cc:54
G4ParticleDefinition * GetParticle(G4int index) const
G4ParticleDefinition * fLastCrossedParticle
bool G4bool
Definition: G4Types.hh:79
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4FastSimulationVector< G4VFastSimulationModel > fInactivatedModels
G4VFastSimulationModel * GetFastSimulationModel(const G4String &modelName, const G4VFastSimulationModel *previousFound, bool &foundPrevious) const
G4VFastSimulationModel * fTriggedFastSimulationModel
void SetCurrentTrack(const G4Track &, const G4Navigator *a=0)
Definition: G4FastTrack.cc:71
static G4TransportationManager * GetTransportationManager()
G4bool InActivateFastSimulationModel(const G4String &)
G4VPhysicalVolume * GetWorldPhysical() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
void SetFastSimulationManager(G4FastSimulationManager *fsm)
Definition: G4Region.cc:128
const G4String & GetParticleName(G4int index) const
G4GLOB_DLL std::ostream G4cout
void ClearFastSimulationManager()
Definition: G4Region.cc:427
G4VPhysicalVolume * GetWorldVolume() const
G4FastSimulationManager(G4Envelope *anEnvelope, G4bool IsUnique=FALSE)
const XML_Char XML_Content * model
Definition: expat.h:151
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
virtual void AtRestDoIt(const G4FastTrack &, G4FastStep &)