Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParallelWorldProcess.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 // $Id: G4ParallelWorldProcess.cc 95193 2016-01-29 08:23:25Z gcosmo $
28 // GEANT4 tag $Name: geant4-09-04-ref-00 $
29 //
30 //
31 
32 #include "G4ios.hh"
35 #include "G4Step.hh"
36 #include "G4StepPoint.hh"
37 #include "G4Navigator.hh"
38 #include "G4VTouchable.hh"
39 #include "G4VPhysicalVolume.hh"
40 #include "G4ParticleChange.hh"
41 #include "G4PathFinder.hh"
43 #include "G4ParticleChange.hh"
44 #include "G4StepPoint.hh"
45 #include "G4FieldTrackUpdator.hh"
46 #include "G4Material.hh"
47 #include "G4ProductionCuts.hh"
48 #include "G4ProductionCutsTable.hh"
49 
50 #include "G4SDManager.hh"
51 #include "G4VSensitiveDetector.hh"
52 
53 G4ThreadLocal G4Step* G4ParallelWorldProcess::fpHyperStep = 0;
54 G4ThreadLocal G4int G4ParallelWorldProcess::nParallelWorlds = 0;
55 G4ThreadLocal G4int G4ParallelWorldProcess::fNavIDHyp = 0;
57 { return fpHyperStep; }
59 { return fNavIDHyp; }
60 
62 G4ParallelWorldProcess(const G4String& processName,G4ProcessType theType)
63 :G4VProcess(processName,theType),fGhostWorld(nullptr),fGhostNavigator(nullptr),
64  fNavigatorID(-1),fFieldTrack('0'),fGhostSafety(0.),fOnBoundary(false),
65  layeredMaterialFlag(false)
66 {
67  SetProcessSubType(491);
68  if(!fpHyperStep) fpHyperStep = new G4Step();
69  iParallelWorld = ++nParallelWorlds;
70 
71  pParticleChange = &aDummyParticleChange;
72 
73  fGhostStep = new G4Step();
74  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
75  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
76 
77  fTransportationManager = G4TransportationManager::GetTransportationManager();
78  fTransportationManager->GetNavigatorForTracking()->SetPushVerbosity(false);
79  fPathFinder = G4PathFinder::GetInstance();
80 
81  fGhostWorldName = "** NotDefined **";
83 
84  if (verboseLevel>0)
85  {
86  G4cout << GetProcessName() << " is created " << G4endl;
87  }
88 }
89 
91 {
92  delete fGhostStep;
93  nParallelWorlds--;
94  if(nParallelWorlds==0)
95  {
96  delete fpHyperStep;
97  fpHyperStep = 0;
98  }
99 }
100 
102 SetParallelWorld(G4String parallelWorldName)
103 {
104  fGhostWorldName = parallelWorldName;
105  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
106  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
107  fGhostNavigator->SetPushVerbosity(false);
108 }
109 
112 {
113  fGhostWorldName = parallelWorld->GetName();
114  fGhostWorld = parallelWorld;
115  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
116  fGhostNavigator->SetPushVerbosity(false);
117 }
118 
120 {
121  if(fGhostNavigator)
122  { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
123  else
124  {
125  G4Exception("G4ParallelWorldProcess::StartTracking",
126  "ProcParaWorld000",FatalException,
127  "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
128  }
129  fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
130 
131  fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
132  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
133  fNewGhostTouchable = fOldGhostTouchable;
134  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
135 
136  fGhostSafety = -1.;
137  fOnBoundary = false;
138  fGhostPreStepPoint->SetStepStatus(fUndefined);
139  fGhostPostStepPoint->SetStepStatus(fUndefined);
140 
141 // G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
142 // if(thePhys)
143 // {
144 // G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
145 // if(ghostMaterial)
146 // { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
147 // }
148 
149  *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
150  if(layeredMaterialFlag)
151  {
152  G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
153  SwitchMaterial(realWorldPostStepPoint);
154  }
155  *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
156 }
157 
158 G4double
160  const G4Track& /*track*/,
162 {
163 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 // At Rest must be registered ONLY for the particle which has other At Rest
165 // process(es).
166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167  *condition = Forced;
168  return DBL_MAX;
169 }
170 
172  const G4Track& track,
173  const G4Step& step)
174 {
175 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
176 // At Rest must be registered ONLY for the particle which has other At Rest
177 // process(es).
178 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
180  G4VSensitiveDetector* aSD = 0;
181  if(fOldGhostTouchable->GetVolume())
182  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
183  fOnBoundary = false;
184  if(aSD)
185  {
186  CopyStep(step);
187  fGhostPreStepPoint->SetSensitiveDetector(aSD);
188 
189  fNewGhostTouchable = fOldGhostTouchable;
190 
191  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
192  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
193  if(fNewGhostTouchable->GetVolume())
194  {
195  fGhostPostStepPoint->SetSensitiveDetector(
196  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
197  }
198  else
199  { fGhostPostStepPoint->SetSensitiveDetector(0); }
200 
201  aSD->Hit(fGhostStep);
202  }
203 
204  pParticleChange->Initialize(track);
205  return pParticleChange;
206 }
207 
208 G4double
210  const G4Track& /*track*/,
211  G4double /*previousStepSize*/,
213 {
214  *condition = StronglyForced;
215  return DBL_MAX;
216 }
217 
219  const G4Track& track,
220  const G4Step& step)
221 {
222  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
223  G4VSensitiveDetector* aSD = 0;
224  if(fOldGhostTouchable->GetVolume())
225  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
226  CopyStep(step);
227  fGhostPreStepPoint->SetSensitiveDetector(aSD);
228 
229  if(fOnBoundary)
230  {
231  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
232  }
233  else
234  {
235  fNewGhostTouchable = fOldGhostTouchable;
236  }
237 
238  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
239  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
240 
241  if(fNewGhostTouchable->GetVolume())
242  {
243  fGhostPostStepPoint->SetSensitiveDetector(
244  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
245  }
246  else
247  { fGhostPostStepPoint->SetSensitiveDetector(0); }
248 
249  G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
250  if(sd)
251  {
252  sd->Hit(fGhostStep);
253  }
254 
255  pParticleChange->Initialize(track);
256  if(layeredMaterialFlag)
257  {
258  G4StepPoint* realWorldPostStepPoint =
259  ((G4Step*)(track.GetStep()))->GetPostStepPoint();
260  SwitchMaterial(realWorldPostStepPoint);
261  }
262  return pParticleChange;
263 }
264 
266  const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
267  G4double& proposedSafety, G4GPILSelection* selection)
268 {
269  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
270  //static ELimited eLimited;
271  ELimited eLimited;
272  ELimited eLim = kUndefLimited;
273 
274  *selection = NotCandidateForSelection;
275  G4double returnedStep = DBL_MAX;
276 
277  if (previousStepSize > 0.)
278  { fGhostSafety -= previousStepSize; }
279  if (fGhostSafety < 0.) fGhostSafety = 0.0;
280 
281  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
282  {
283  // I have no chance to limit
284  returnedStep = currentMinimumStep;
285  fOnBoundary = false;
286  proposedSafety = fGhostSafety - currentMinimumStep;
287  eLim = kDoNot;
288  }
289  else
290  {
291  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
292 
293 #ifdef G4DEBUG_PARALLEL_WORLD_PROCESS
294  if( verboseLevel > 0 ){
295  int localVerb = verboseLevel-1;
296 
297  if( localVerb == 1 ) {
298  G4cout << " Pll Wrl proc::AlongStepGPIL " << this->GetProcessName() << G4endl;
299  }else if( localVerb > 1 ) {
300  G4cout << "----------------------------------------------" << G4endl;
301  G4cout << " ParallelWorldProcess: field Track set to : " << G4endl;
302  G4cout << "----------------------------------------------" << G4endl;
303  G4cout << fFieldTrack << G4endl;
304  G4cout << "----------------------------------------------" << G4endl;
305  }
306  }
307 #endif
308 
309  returnedStep
310  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
311  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
312  endTrack,track.GetVolume());
313  if(eLimited == kDoNot)
314  {
315  fOnBoundary = false;
316  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
317  }
318  else
319  {
320  fOnBoundary = true;
321  // fGhostSafetyEnd = 0.0; // At end-point of expected step only
322  }
323  proposedSafety = fGhostSafety;
324  if(eLimited == kUnique || eLimited == kSharedOther) {
325  *selection = CandidateForSelection;
326  }
327  else if (eLimited == kSharedTransport) {
328  returnedStep *= (1.0 + 1.0e-9);
329  }
330  eLim = eLimited;
331  }
332 
333  if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
334  if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
335  return returnedStep;
336 }
337 
339  const G4Track& track, const G4Step& )
340 {
341  pParticleChange->Initialize(track);
342  return pParticleChange;
343 }
344 
345 void G4ParallelWorldProcess::CopyStep(const G4Step & step)
346 {
347  G4StepStatus prevStat = fGhostPostStepPoint->GetStepStatus();
348 
349  fGhostStep->SetTrack(step.GetTrack());
350  fGhostStep->SetStepLength(step.GetStepLength());
351  fGhostStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
353  fGhostStep->SetControlFlag(step.GetControlFlag());
354 
355  *fGhostPreStepPoint = *(step.GetPreStepPoint());
356  *fGhostPostStepPoint = *(step.GetPostStepPoint());
357 
358  fGhostPreStepPoint->SetStepStatus(prevStat);
359  if(fOnBoundary)
360  { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
361  else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
362  { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
363 
364  if(iParallelWorld==1)
365  {
366  G4StepStatus prevStatHyp = fpHyperStep->GetPostStepPoint()->GetStepStatus();
367 
368  fpHyperStep->SetTrack(step.GetTrack());
369  fpHyperStep->SetStepLength(step.GetStepLength());
370  fpHyperStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
372  fpHyperStep->SetControlFlag(step.GetControlFlag());
373 
374  *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
375  *(fpHyperStep->GetPostStepPoint()) = *(step.GetPostStepPoint());
376 
377  fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
378  }
379 
380  if(fOnBoundary)
381  { fpHyperStep->GetPostStepPoint()->SetStepStatus(fGeomBoundary); }
382 }
383 
384 void G4ParallelWorldProcess::SwitchMaterial(G4StepPoint* realWorldStepPoint)
385 {
386  if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
387  G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
388  if(thePhys)
389  {
390  G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
391  if(ghostMaterial)
392  {
393  G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
394  G4ProductionCuts* prodCuts =
395  realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
396  if(ghostRegion)
397  {
398  G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
399  if(ghostProdCuts) prodCuts = ghostProdCuts;
400  }
401  const G4MaterialCutsCouple* ghostMCCouple =
403  ->GetMaterialCutsCouple(ghostMaterial,prodCuts);
404  if(ghostMCCouple)
405  {
406  realWorldStepPoint->SetMaterial(ghostMaterial);
407  realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
408  *(fpHyperStep->GetPostStepPoint()) = *(fGhostPostStepPoint);
409  fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
410  fpHyperStep->GetPostStepPoint()->SetMaterialCutsCouple(ghostMCCouple);
411  }
412  else
413  {
414  G4cout << "!!! MaterialCutsCouple is not found for "
415  << ghostMaterial->GetName() << "." << G4endl
416  << " Material in real world ("
417  << realWorldStepPoint->GetMaterial()->GetName()
418  << ") is used." << G4endl;
419  }
420  }
421  }
422 }
423 
425 {
426  G4int pdgCode = partDef->GetPDGEncoding();
427  if(pdgCode==0)
428  {
429  G4String partName = partDef->GetParticleName();
430  if(partName=="opticalphoton") return false;
431  if(partName=="geantino") return false;
432  if(partName=="chargedgeantino") return false;
433  }
434  else
435  {
436  if(pdgCode==22) return false; // gamma
437  if(pdgCode==11) return false; // electron
438  if(pdgCode==2212) return false; // proton
439  if(pdgCode==-12) return false; // anti_nu_e
440  if(pdgCode==12) return false; // nu_e
441  if(pdgCode==-14) return false; // anti_nu_mu
442  if(pdgCode==14) return false; // nu_mu
443  if(pdgCode==-16) return false; // anti_nu_tau
444  if(pdgCode==16) return false; // nu_tau
445  }
446  return true;
447 }
448 
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4double condition(const G4ErrorSymMatrix &m)
void SetStepLength(G4double value)
virtual void Initialize(const G4Track &)
G4ProductionCuts * GetProductionCuts() const
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetTrack(G4Track *value)
G4Material * GetMaterial() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetStepLength() const
G4bool IsAtRestRequired(G4ParticleDefinition *)
G4StepStatus GetStepStatus() const
const G4String & GetName() const
Definition: G4Material.hh:178
G4Material * GetMaterial() const
ELimited
G4double GetNonIonizingEnergyDeposit() const
const G4ThreeVector & GetPosition() const
G4Navigator * GetNavigatorForTracking() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
G4SteppingControl GetControlFlag() const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
G4Region * GetRegion() const
#define G4ThreadLocal
Definition: tls.hh:89
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetStepStatus(const G4StepStatus aValue)
const G4String & GetParticleName() const
G4StepStatus
Definition: G4StepStatus.hh:49
static G4ParallelWorldProcessStore * GetInstance()
G4StepPoint * GetPreStepPoint() const
void SetParallelWorld(G4String parallelWorldName)
void SetSensitiveDetector(G4VSensitiveDetector *)
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
const G4String & GetName() const
void SetControlFlag(G4SteppingControl StepControlFlag)
bool G4bool
Definition: G4Types.hh:79
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void SetParallelWorld(G4ParallelWorldProcess *proc, G4String parallelWorldName)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4bool Hit(G4Step *aStep)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTotalEnergyDeposit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
static G4ProductionCutsTable * GetProductionCutsTable()
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4int ActivateNavigator(G4Navigator *aNavigator)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4Navigator * GetNavigator(const G4String &worldName)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
void SetNonIonizingEnergyDeposit(G4double value)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VPhysicalVolume * GetVolume() const
void SetPushVerbosity(G4bool mode)
#define G4endl
Definition: G4ios.hh:61
void SetTotalEnergyDeposit(G4double value)
G4VSensitiveDetector * GetSensitiveDetector() const
void SetMaterial(G4Material *)
G4ParallelWorldProcess(const G4String &processName="ParaWorld", G4ProcessType theType=fParallel)
double G4double
Definition: G4Types.hh:76
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4ForceCondition
G4Track * GetTrack() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const
static void Update(G4FieldTrack *, const G4Track *)
static const G4Step * GetHyperStep()
G4GPILSelection
G4ProcessType