source: trunk/source/processes/parameterisation/src/G4FastSimulationManager.cc @ 1337

Last change on this file since 1337 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 13.5 KB
Line 
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,v 1.13 2007/05/11 13:50:20 mverderi Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30//---------------------------------------------------------------
31//
32//  G4FastSimulationManager.cc
33//
34//  Description:
35//    Manages the Fast Simulation models attached to a envelope.
36//
37//  History:
38//    Oct 97: Verderi && MoraDeFreitas - First Implementation.
39//    ...
40//    May 07: Move to parallel world scheme
41//
42//---------------------------------------------------------------
43
44#include "G4FastSimulationManager.hh"
45#include "G4GlobalFastSimulationManager.hh"
46#include "G4PVPlacement.hh"
47#include "G4TransportationManager.hh"
48
49// --------------------------------------------------
50// Constructor with envelope and IsUnique flag :
51// --------------------------------------------------
52//
53G4FastSimulationManager::
54G4FastSimulationManager(G4Envelope *anEnvelope,
55                        G4bool IsUnique) :
56  fFastTrack(anEnvelope,IsUnique),fTriggedFastSimulationModel(0),
57  fLastCrossedParticle(0)
58{
59  // Communicates to the region that it becomes a
60  // envelope and with this fast simulation manager.
61  anEnvelope->SetFastSimulationManager(this);
62
63  // Add itself to the GlobalFastSimulationManager
64  G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->
65    AddFastSimulationManager(this);
66}
67
68// -----------
69// Destructor:
70// -----------
71G4FastSimulationManager::~G4FastSimulationManager()
72{
73  //
74  // Check out the Envelope about this pointer. If in use,
75  // resets the Logical Volume IsEnvelope flag to avoid clash.
76  //
77  if(fFastTrack.GetEnvelope()->GetFastSimulationManager()==this)
78    fFastTrack.GetEnvelope()->ClearFastSimulationManager();
79  // Remove itself from the GlobalFastSimulationManager
80  G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->
81    RemoveFastSimulationManager(this);
82}
83
84// ---------------------------------------
85// Methods to activate/inactivate models
86//----------------------------------------
87
88G4bool
89G4FastSimulationManager::ActivateFastSimulationModel(const G4String& aName) 
90{
91  size_t iModel;
92
93  // If the model is already active, do nothing.
94  for (iModel=0; iModel<ModelList.size(); iModel++)
95    if(ModelList[iModel]->GetName() == aName)
96      return true;
97 
98  // Look for in the fInactivatedModels list, if found push_back it back to
99  // the ModelList
100  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
101    if(fInactivatedModels[iModel]->GetName() == aName) {
102      ModelList.
103        push_back (fInactivatedModels.removeAt(iModel));
104      // forces the fApplicableModelList to be rebuild
105      fLastCrossedParticle=0;
106      return true;
107    }
108  return false;
109}
110
111G4bool
112G4FastSimulationManager::InActivateFastSimulationModel(const G4String& aName)
113{
114  // Look for in the ModelList, if found remove from it and keep the pointer
115  // on the fInactivatedModels list.
116  for (size_t iModel=0; iModel<ModelList.size(); iModel++)
117    if(ModelList[iModel]->GetName() == aName) {
118      fInactivatedModels.
119        push_back (ModelList.removeAt(iModel));
120      // forces the fApplicableModelList to be rebuild
121      fLastCrossedParticle=0;
122      return true;
123    }
124  return false;
125}
126
127G4VFastSimulationModel* 
128G4FastSimulationManager::GetFastSimulationModel(const G4String& modelName,
129                                                const G4VFastSimulationModel* previousFound,
130                                                bool &foundPrevious) const
131{
132  G4VFastSimulationModel* model = 0;
133  for (size_t iModel=0; iModel<ModelList.size(); iModel++)
134    {
135      if(ModelList[iModel]->GetName() == modelName)
136        {
137          if (previousFound == 0)
138            {
139              model = ModelList[iModel];
140              break;
141            }
142          else
143            {
144              if (ModelList[iModel] == previousFound)
145                {
146                  foundPrevious = true;
147                  continue;
148                }
149              if (foundPrevious)
150                {
151                  model = ModelList[iModel];
152                  break;
153                }
154            }
155        }
156    }
157  return model;
158}
159
160
161//------------------------------------------------------------------
162// Interface trigger method for the G4ParameterisationManagerProcess
163//------------------------------------------------------------------
164//   G4bool GetFastSimulationManagerTrigger(const G4Track &);
165//
166//    This method is used to interface the G4FastSimulationManagerProcess
167//    with the user Fast Simulation Models. It's called when the particle
168//    is inside the envelope.
169//
170//    It :
171//
172//      1) initialises the private members (fFastTrack and so
173//         on);
174//      2) loops on the IsApplicable() methods to find out the
175//         ones should be applied.
176//      2) for these, loops on the ModelTrigger() methods to find out
177//         perhaps one that must be applied just now.
178//
179//    If the a Fast Simulation Model is triggered then it returns
180//    true, false otherwise.
181//
182//-----------------------------------------------------------
183G4bool
184G4FastSimulationManager::
185PostStepGetFastSimulationManagerTrigger(const G4Track& track,
186                                        const G4Navigator* theNavigator)
187{
188  size_t iModel;
189 
190  // If particle type changed re-build the fApplicableModelList.
191  if(fLastCrossedParticle!=track.GetDefinition()) {
192    fLastCrossedParticle=track.GetDefinition();
193    fApplicableModelList.clear();
194    // If Model List is empty, do nothing !
195    if(ModelList.size()==0) return false;
196    for (iModel=0; iModel<ModelList.size(); iModel++)
197      if(ModelList[iModel]->IsApplicable(*(track.GetDefinition())))
198        fApplicableModelList.push_back (ModelList[iModel]);
199  }
200
201  // If Applicable Model List is empty, do nothing !
202  if(fApplicableModelList.size()==0) return false;
203
204  // -- Register current track
205  fFastTrack.SetCurrentTrack(track,theNavigator);
206
207  // tests if particle are on the boundary and leaving,
208  // in this case do nothing !
209  if(fFastTrack.OnTheBoundaryButExiting()) return false;
210 
211  // Loops on the ModelTrigger() methods
212  for (iModel=0; iModel<fApplicableModelList.size(); iModel++)
213   
214    //---------------------------------------------------
215    // Asks the ModelTrigger method if it must be trigged now.
216    //---------------------------------------------------
217   
218    if(fApplicableModelList[iModel]->ModelTrigger(fFastTrack)) {
219      //--------------------------------------------------
220      // The model will be applied. Initializes the G4FastStep
221      // with the current state of the G4Track and
222      // same usefull parameters.
223      // In particular it does SetLocalEnergyDeposit(0.0).
224      //--------------------------------------------------     
225      fFastStep.Initialize(fFastTrack);
226     
227      // Keeps the FastSimulationModel pointer to call the
228      // DoIt() method.
229      fTriggedFastSimulationModel=fApplicableModelList[iModel];
230      return true;
231    }
232
233  //--------------------------------------------
234  // Nobody asks to gain control, returns false
235  //--------------------------------------------
236  return false;
237}
238
239G4VParticleChange* G4FastSimulationManager::InvokePostStepDoIt() 
240{
241  //  const G4FastTrack& parFastTrack=fFastTrack;
242  fTriggedFastSimulationModel->DoIt(fFastTrack,fFastStep);
243  return &fFastStep;
244}
245
246// -------------------------------------------------------------
247// -- Mostly the same as above, in the case of AtRest particles:
248// -------------------------------------------------------------
249G4bool
250G4FastSimulationManager::AtRestGetFastSimulationManagerTrigger(const G4Track& track,
251                                                               const G4Navigator* theNavigator)
252{
253  size_t iModel;
254 
255  // If particle type changed re-build the fApplicableModelList.
256  if(fLastCrossedParticle!=track.GetDefinition()) {
257    fLastCrossedParticle=track.GetDefinition();
258    fApplicableModelList.clear();
259    // If Model List is empty, do nothing !
260    if(ModelList.size()==0) return false;
261    for (iModel=0; iModel<ModelList.size(); iModel++)
262      if(ModelList[iModel]->IsApplicable(*(track.GetDefinition())))
263        fApplicableModelList.push_back (ModelList[iModel]);
264  }
265 
266  // If Applicable Model List is empty, do nothing !
267  if(fApplicableModelList.size()==0) return false;
268
269  // -- Register current track
270  fFastTrack.SetCurrentTrack(track,theNavigator);
271 
272  // -- (note: compared to the PostStepGetFastSimulationManagerTrigger,
273  // --  the test to see if the particle is on the boundary but leaving
274  // --  is irrelevant here)
275 
276  // Loops on the models to see if one of them wants to trigger:
277  for (iModel=0; iModel < fApplicableModelList.size(); iModel++)
278    if(fApplicableModelList[iModel]->AtRestModelTrigger(fFastTrack))
279      {
280        fFastStep.Initialize(fFastTrack);
281        fTriggedFastSimulationModel=fApplicableModelList[iModel];
282        return true;
283      }
284 
285  //--------------------------------------------
286  // Nobody asks to gain control, returns false
287  //--------------------------------------------
288  return false;
289}
290
291G4VParticleChange* G4FastSimulationManager::InvokeAtRestDoIt() 
292{
293  fTriggedFastSimulationModel->AtRestDoIt(fFastTrack,fFastStep);
294  return &fFastStep;
295}
296
297void 
298G4FastSimulationManager::ListTitle() const
299{
300  G4cout << fFastTrack.GetEnvelope()->GetName();
301  //  if(GhostPlacements.size()!=0) G4cout << " (ghost)";
302  if (fFastTrack.GetEnvelope()->GetWorldPhysical() == G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume()) G4cout << " (mass geom.)";
303  else G4cout << " (// geom.)";
304                                                                                                                                                         
305}
306
307void 
308G4FastSimulationManager::ListModels() const
309{
310  size_t iModel;
311
312  G4cout << "Current Models for the ";
313  ListTitle();
314  G4cout << " envelope:\n";
315
316  for (iModel=0; iModel<ModelList.size(); iModel++) 
317    G4cout << "   " << ModelList[iModel]->GetName() << "\n";
318
319  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
320    G4cout << "   " << fInactivatedModels[iModel]->GetName() 
321           << "(inactivated)\n";
322}
323
324void 
325G4FastSimulationManager::ListModels(const G4String& aName) const
326{
327  size_t iModel;
328  G4int titled = 0;
329  G4ParticleTable* theParticleTable=
330    G4ParticleTable::GetParticleTable();
331 
332  // Active Models
333  for (iModel=0; iModel<ModelList.size(); iModel++)
334    if(ModelList[iModel]->GetName() == aName ||
335       aName == "all" ) {
336      if(!(titled++)){
337        G4cout << "In the envelope ";
338        ListTitle();
339        G4cout << ",\n";
340      }
341      G4cout << "  the model " << ModelList[iModel]->GetName()
342             << " is applicable for :\n     ";
343     
344      G4int list_started=0;
345      for (G4int iParticle=0; iParticle<theParticleTable->entries(); 
346           iParticle++)
347        if(ModelList[iModel]->
348           IsApplicable(*(theParticleTable->
349                          GetParticle(iParticle)))) {
350          if(list_started++) G4cout << ", ";
351          G4cout << theParticleTable->
352            GetParticle(iParticle)->GetParticleName();
353        }
354      G4cout <<G4endl;
355    }
356 
357  // Inactive Models
358  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
359    if(fInactivatedModels[iModel]->GetName() == aName ||
360       aName == "all" ) {
361      if(!(titled++)){
362        G4cout << "In the envelope ";
363        ListTitle();
364        G4cout << ",\n";
365      }
366      G4cout << "  the model " << fInactivatedModels[iModel]->GetName()
367             << " (inactivated) is applicable for :\n     ";
368     
369      G4int list_started=0;
370      for (G4int iParticle=0; iParticle<theParticleTable->entries(); 
371           iParticle++)
372        if(fInactivatedModels[iModel]->
373           IsApplicable(*(theParticleTable->
374                          GetParticle(iParticle)))) {
375          if(list_started++) G4cout << ", ";
376          G4cout << theParticleTable->
377            GetParticle(iParticle)->GetParticleName();
378        }
379      G4cout <<G4endl;
380    }
381}
382
383void 
384G4FastSimulationManager::ListModels(const G4ParticleDefinition* aPD) const
385{
386  size_t iModel;
387  G4bool unique=true;
388 
389  // Active Models
390  for (iModel=0; iModel<ModelList.size(); iModel++)
391    if(ModelList[iModel]->IsApplicable(*aPD)) {
392      G4cout << "Envelope ";
393      ListTitle();
394      G4cout << ", Model " 
395             << ModelList[iModel]->GetName() 
396             << "." << G4endl;
397    }
398  // inactive Models
399  for (iModel=0; iModel<fInactivatedModels.size(); iModel++)
400    if(fInactivatedModels[iModel]->IsApplicable(*aPD)) {
401      G4cout << "Envelope ";
402      ListTitle();
403      G4cout << ", Model " 
404             << fInactivatedModels[iModel]->GetName() 
405             << " (inactivated)." << G4endl;
406    }
407 
408  if(!unique)
409    G4cout << "\a\n >>>>>>Warning: two or more Models for the same "
410           << "particle type attached to the same envelope!"
411           << G4endl;
412  unique=false;
413}
Note: See TracBrowser for help on using the repository browser.