Geant4  10.00.p02
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 76935 2013-11-19 09:46:31Z gcosmo $
28 // GEANT4 tag $Name: geant4-09-04-ref-00 $
29 //
30 //
31 
32 #include "G4ios.hh"
34 #include "G4Step.hh"
35 #include "G4StepPoint.hh"
36 #include "G4Navigator.hh"
37 #include "G4VTouchable.hh"
38 #include "G4VPhysicalVolume.hh"
39 #include "G4ParticleChange.hh"
40 #include "G4PathFinder.hh"
42 #include "G4ParticleChange.hh"
43 #include "G4StepPoint.hh"
44 #include "G4FieldTrackUpdator.hh"
45 #include "G4Material.hh"
46 #include "G4ProductionCuts.hh"
47 #include "G4ProductionCutsTable.hh"
48 
49 #include "G4SDManager.hh"
50 #include "G4VSensitiveDetector.hh"
51 
56 { return fpHyperStep; }
58 { return fNavIDHyp; }
59 
61 G4ParallelWorldProcess(const G4String& processName,G4ProcessType theType)
62 :G4VProcess(processName,theType), fGhostNavigator(0), fNavigatorID(-1),
63  fFieldTrack('0'),layeredMaterialFlag(false)
64 {
65  if(!fpHyperStep) fpHyperStep = new G4Step();
67 
69 
70  fGhostStep = new G4Step();
73 
77 
78  if (verboseLevel>0)
79  {
80  G4cout << GetProcessName() << " is created " << G4endl;
81  }
82 }
83 
85 {
86  delete fGhostStep;
88  if(nParallelWorlds==0)
89  {
90  delete fpHyperStep;
91  fpHyperStep = 0;
92  }
93 }
94 
96 SetParallelWorld(G4String parallelWorldName)
97 {
98  fGhostWorldName = parallelWorldName;
102 }
103 
106 {
107  fGhostWorldName = parallelWorld->GetName();
108  fGhostWorld = parallelWorld;
111 }
112 
114 {
115  if(fGhostNavigator)
117  else
118  {
119  G4Exception("G4ParallelWorldProcess::StartTracking",
120  "ProcParaWorld000",FatalException,
121  "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
122  }
124 
129 
130  fGhostSafety = -1.;
131  fOnBoundary = false;
134 
135 // G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
136 // if(thePhys)
137 // {
138 // G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
139 // if(ghostMaterial)
140 // { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
141 // }
142 
145  {
146  G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
147  SwitchMaterial(realWorldPostStepPoint);
148  }
150 }
151 
152 G4double
154  const G4Track& /*track*/,
156 {
157 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 // At Rest must be registered ONLY for the particle which has other At Rest
159 // process(es).
160 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
161  *condition = Forced;
162  return DBL_MAX;
163 }
164 
166  const G4Track& track,
167  const G4Step& step)
168 {
169 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
170 // At Rest must be registered ONLY for the particle which has other At Rest
171 // process(es).
172 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
174  G4VSensitiveDetector* aSD = 0;
177  fOnBoundary = false;
178  if(aSD)
179  {
180  CopyStep(step);
182 
184 
188  {
191  }
192  else
194 
195  aSD->Hit(fGhostStep);
196  }
197 
198  pParticleChange->Initialize(track);
199  return pParticleChange;
200 }
201 
202 G4double
204  const G4Track& /*track*/,
205  G4double /*previousStepSize*/,
207 {
208  *condition = StronglyForced;
209  return DBL_MAX;
210 }
211 
213  const G4Track& track,
214  const G4Step& step)
215 {
217  G4VSensitiveDetector* aSD = 0;
220  CopyStep(step);
222 
223  if(fOnBoundary)
224  {
226  }
227  else
228  {
230  }
231 
234 
236  {
239  }
240  else
242 
244  if(sd)
245  {
246  sd->Hit(fGhostStep);
247  }
248 
249  pParticleChange->Initialize(track);
251  {
252  G4StepPoint* realWorldPostStepPoint =
253  ((G4Step*)(track.GetStep()))->GetPostStepPoint();
254  SwitchMaterial(realWorldPostStepPoint);
255  }
256  return pParticleChange;
257 }
258 
260  const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
261  G4double& proposedSafety, G4GPILSelection* selection)
262 {
263  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
264  //static ELimited eLimited;
265  ELimited eLimited;
266  ELimited eLim = kUndefLimited;
267 
268  *selection = NotCandidateForSelection;
269  G4double returnedStep = DBL_MAX;
270 
271  if (previousStepSize > 0.)
272  { fGhostSafety -= previousStepSize; }
273  if (fGhostSafety < 0.) fGhostSafety = 0.0;
274 
275  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
276  {
277  // I have no chance to limit
278  returnedStep = currentMinimumStep;
279  fOnBoundary = false;
280  proposedSafety = fGhostSafety - currentMinimumStep;
281  eLim = kDoNot;
282  }
283  else
284  {
286  returnedStep
287  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
288  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
289  endTrack,track.GetVolume());
290  if(eLimited == kDoNot)
291  {
292  fOnBoundary = false;
293  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
294  }
295  else
296  {
297  fOnBoundary = true;
298  }
299  proposedSafety = fGhostSafety;
300  if(eLimited == kUnique || eLimited == kSharedOther) {
301  *selection = CandidateForSelection;
302  }
303  else if (eLimited == kSharedTransport) {
304  returnedStep *= (1.0 + 1.0e-9);
305  }
306  eLim = eLimited;
307  }
308 
310  if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
311  return returnedStep;
312 }
313 
315  const G4Track& track, const G4Step& )
316 {
317  pParticleChange->Initialize(track);
318  return pParticleChange;
319 }
320 
322 {
324 
325  fGhostStep->SetTrack(step.GetTrack());
330 
331  *fGhostPreStepPoint = *(step.GetPreStepPoint());
333 
335  if(fOnBoundary)
339 
340  if(iParallelWorld==1)
341  {
343 
344  fpHyperStep->SetTrack(step.GetTrack());
349 
352 
353  fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
354  }
355 
356  if(fOnBoundary)
358 }
359 
361 {
362  if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
364  if(thePhys)
365  {
366  G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
367  if(ghostMaterial)
368  {
369  G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
370  G4ProductionCuts* prodCuts =
371  realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
372  if(ghostRegion)
373  {
374  G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
375  if(ghostProdCuts) prodCuts = ghostProdCuts;
376  }
377  const G4MaterialCutsCouple* ghostMCCouple =
379  ->GetMaterialCutsCouple(ghostMaterial,prodCuts);
380  if(ghostMCCouple)
381  {
382  realWorldStepPoint->SetMaterial(ghostMaterial);
383  realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
385  fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
387  }
388  else
389  {
390  G4cout << "!!! MaterialCutsCouple is not found for "
391  << ghostMaterial->GetName() << "." << G4endl
392  << " Material in real world ("
393  << realWorldStepPoint->GetMaterial()->GetName()
394  << ") is used." << G4endl;
395  }
396  }
397  }
398 }
399 
401 {
402  G4int pdgCode = partDef->GetPDGEncoding();
403  if(pdgCode==0)
404  {
405  G4String partName = partDef->GetParticleName();
406  if(partName=="opticalphoton") return false;
407  if(partName=="geantino") return false;
408  if(partName=="chargedgeantino") return false;
409  }
410  else
411  {
412  if(pdgCode==22) return false; // gamma
413  if(pdgCode==11) return false; // electron
414  if(pdgCode==2212) return false; // proton
415  if(pdgCode==-12) return false; // anti_nu_e
416  if(pdgCode==12) return false; // nu_e
417  if(pdgCode==-14) return false; // anti_nu_mu
418  if(pdgCode==14) return false; // nu_mu
419  if(pdgCode==-16) return false; // anti_nu_tau
420  if(pdgCode==16) return false; // nu_tau
421  }
422  return true;
423 }
424 
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
void CopyStep(const G4Step &step)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetTrack(G4Track *value)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetStepLength() const
G4Material * GetMaterial() const
G4bool IsAtRestRequired(G4ParticleDefinition *)
G4StepStatus GetStepStatus() const
const G4String & GetName() const
Definition: G4Material.hh:176
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
G4VPhysicalVolume * fGhostWorld
#define G4ThreadLocal
Definition: tls.hh:52
G4TouchableHandle fNewGhostTouchable
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
G4StepPoint * GetPreStepPoint() const
void SetParallelWorld(G4String parallelWorldName)
void SetSensitiveDetector(G4VSensitiveDetector *)
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
static G4ThreadLocal G4int nParallelWorlds
const G4String & GetName() const
void SetControlFlag(G4SteppingControl StepControlFlag)
bool G4bool
Definition: G4Types.hh:79
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
static G4ThreadLocal G4int fNavIDHyp
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4bool Hit(G4Step *aStep)
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()
void SwitchMaterial(G4StepPoint *)
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
G4VParticleChange aDummyParticleChange
G4StepPoint * GetPostStepPoint() const
void SetNonIonizingEnergyDeposit(G4double value)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VPhysicalVolume * GetVolume() const
G4ParallelWorldProcess(const G4String &processName="ParaWorld", G4ProcessType theType=fParameterisation)
void SetPushVerbosity(G4bool mode)
#define G4endl
Definition: G4ios.hh:61
void SetTotalEnergyDeposit(G4double value)
G4VSensitiveDetector * GetSensitiveDetector() const
void SetMaterial(G4Material *)
double G4double
Definition: G4Types.hh:76
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
G4TouchableHandle fOldGhostTouchable
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4ForceCondition
G4Track * GetTrack() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
static G4ThreadLocal G4Step * fpHyperStep
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
G4TransportationManager * fTransportationManager
const G4TouchableHandle & GetTouchableHandle() const
static void Update(G4FieldTrack *, const G4Track *)
static const G4Step * GetHyperStep()
G4GPILSelection
G4ProcessType