Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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"
45#include "G4Material.hh"
46#include "G4ProductionCuts.hh"
48
49#include "G4SDManager.hh"
51
52G4Step* G4ParallelWorldProcess::fpHyperStep = 0;
53G4int G4ParallelWorldProcess::nParallelWorlds = 0;
54G4int G4ParallelWorldProcess::fNavIDHyp = 0;
56{ return fpHyperStep; }
58{ return fNavIDHyp; }
59
61G4ParallelWorldProcess(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();
66 iParallelWorld = ++nParallelWorlds;
67
68 pParticleChange = &aDummyParticleChange;
69
70 fGhostStep = new G4Step();
71 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
72 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
73
75 fPathFinder = G4PathFinder::GetInstance();
76
77 if (verboseLevel>0)
78 {
79 G4cout << GetProcessName() << " is created " << G4endl;
80 }
81}
82
84{
85 delete fGhostStep;
86 nParallelWorlds--;
87 if(nParallelWorlds==0)
88 {
89 delete fpHyperStep;
90 fpHyperStep = 0;
91 }
92}
93
95SetParallelWorld(G4String parallelWorldName)
96{
97 fGhostWorldName = parallelWorldName;
98 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
99 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
100}
101
104{
105 fGhostWorldName = parallelWorld->GetName();
106 fGhostWorld = parallelWorld;
107 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
108}
109
111{
112 if(fGhostNavigator)
113 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
114 else
115 {
116 G4Exception("G4ParallelWorldProcess::StartTracking",
117 "ProcParaWorld000",FatalException,
118 "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
119 }
120 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
121
122 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
123 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
124 fNewGhostTouchable = fOldGhostTouchable;
125 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
126
127 fGhostSafety = -1.;
128 fOnBoundary = false;
129 fGhostPreStepPoint->SetStepStatus(fUndefined);
130 fGhostPostStepPoint->SetStepStatus(fUndefined);
131
132// G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
133// G4cout << "======= G4ParallelWorldProcess::StartTracking() =======" << G4endl;
134// if(thePhys)
135// {
136// G4cout << " --- Parallel world volume : " << thePhys->GetName() << G4endl;
137// G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
138// if(ghostMaterial)
139// { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
140// }
141
142 *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
143 if(layeredMaterialFlag)
144 {
145 G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
146 SwitchMaterial(realWorldPostStepPoint);
147 }
148 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
149}
150
153 const G4Track& /*track*/,
155{
156//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
157// At Rest must be registered ONLY for the particle which has other At Rest
158// process(es).
159//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160 *condition = Forced;
161 return DBL_MAX;
162}
163
165 const G4Track& track,
166 const G4Step& step)
167{
168//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169// At Rest must be registered ONLY for the particle which has other At Rest
170// process(es).
171//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
172 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
173 G4VSensitiveDetector* aSD = 0;
174 if(fOldGhostTouchable->GetVolume())
175 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
176 fOnBoundary = false;
177 if(aSD)
178 {
179 CopyStep(step);
180 fGhostPreStepPoint->SetSensitiveDetector(aSD);
181
182 fNewGhostTouchable = fOldGhostTouchable;
183
184 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
185 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
186 if(fNewGhostTouchable->GetVolume())
187 {
188 fGhostPostStepPoint->SetSensitiveDetector(
189 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
190 }
191 else
192 { fGhostPostStepPoint->SetSensitiveDetector(0); }
193
194 aSD->Hit(fGhostStep);
195 }
196
198 return pParticleChange;
199}
200
203 const G4Track& /*track*/,
204 G4double /*previousStepSize*/,
206{
208 return DBL_MAX;
209}
210
212 const G4Track& track,
213 const G4Step& step)
214{
215 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
216 G4VSensitiveDetector* aSD = 0;
217 if(fOldGhostTouchable->GetVolume())
218 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
219 CopyStep(step);
220 fGhostPreStepPoint->SetSensitiveDetector(aSD);
221
222 if(fOnBoundary)
223 {
224 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
225 }
226 else
227 {
228 fNewGhostTouchable = fOldGhostTouchable;
229 }
230
231 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
232 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
233
234 if(fNewGhostTouchable->GetVolume())
235 {
236 fGhostPostStepPoint->SetSensitiveDetector(
237 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
238 }
239 else
240 { fGhostPostStepPoint->SetSensitiveDetector(0); }
241
242 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
243 if(sd)
244 {
245 sd->Hit(fGhostStep);
246 }
247
249 if(layeredMaterialFlag)
250 {
251 G4StepPoint* realWorldPostStepPoint =
252 ((G4Step*)(track.GetStep()))->GetPostStepPoint();
253 SwitchMaterial(realWorldPostStepPoint);
254 }
255 return pParticleChange;
256}
257
259 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
260 G4double& proposedSafety, G4GPILSelection* selection)
261{
262 static G4FieldTrack endTrack('0');
263 //static ELimited eLimited;
264 ELimited eLimited;
265 ELimited eLim = kUndefLimited;
266
267 *selection = NotCandidateForSelection;
268 G4double returnedStep = DBL_MAX;
269
270 if (previousStepSize > 0.)
271 { fGhostSafety -= previousStepSize; }
272 if (fGhostSafety < 0.) fGhostSafety = 0.0;
273
274 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
275 {
276 // I have no chance to limit
277 returnedStep = currentMinimumStep;
278 fOnBoundary = false;
279 proposedSafety = fGhostSafety - currentMinimumStep;
280 eLim = kDoNot;
281 }
282 else
283 {
284 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
285 returnedStep
286 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
287 track.GetCurrentStepNumber(),fGhostSafety,eLimited,
288 endTrack,track.GetVolume());
289 if(eLimited == kDoNot)
290 {
291 fOnBoundary = false;
292 fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
293 }
294 else
295 {
296 fOnBoundary = true;
297 }
298 proposedSafety = fGhostSafety;
299 if(eLimited == kUnique || eLimited == kSharedOther) {
300 *selection = CandidateForSelection;
301 }
302 else if (eLimited == kSharedTransport) {
303 returnedStep *= (1.0 + 1.0e-9);
304 }
305 eLim = eLimited;
306 }
307
308 if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
309 if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
310 return returnedStep;
311}
312
314 const G4Track& track, const G4Step& )
315{
317 return pParticleChange;
318}
319
320void G4ParallelWorldProcess::CopyStep(const G4Step & step)
321{
322 G4StepStatus prevStat = fGhostPostStepPoint->GetStepStatus();
323
324 fGhostStep->SetTrack(step.GetTrack());
325 fGhostStep->SetStepLength(step.GetStepLength());
328 fGhostStep->SetControlFlag(step.GetControlFlag());
329
330 *fGhostPreStepPoint = *(step.GetPreStepPoint());
331 *fGhostPostStepPoint = *(step.GetPostStepPoint());
332
333 fGhostPreStepPoint->SetStepStatus(prevStat);
334 if(fOnBoundary)
335 { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
336 else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
337 { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
338
339 if(iParallelWorld==1)
340 {
341 G4StepStatus prevStatHyp = fpHyperStep->GetPostStepPoint()->GetStepStatus();
342
343 fpHyperStep->SetTrack(step.GetTrack());
344 fpHyperStep->SetStepLength(step.GetStepLength());
345 fpHyperStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
347 fpHyperStep->SetControlFlag(step.GetControlFlag());
348
349 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
350 *(fpHyperStep->GetPostStepPoint()) = *(step.GetPostStepPoint());
351
352 fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
353 }
354
355 if(fOnBoundary)
356 { fpHyperStep->GetPostStepPoint()->SetStepStatus(fGeomBoundary); }
357}
358
359void G4ParallelWorldProcess::SwitchMaterial(G4StepPoint* realWorldStepPoint)
360{
361 if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
362 G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
363 if(thePhys)
364 {
365 G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
366 if(ghostMaterial)
367 {
368 G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
369 G4ProductionCuts* prodCuts =
370 realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
371 if(ghostRegion)
372 {
373 G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
374 if(ghostProdCuts) prodCuts = ghostProdCuts;
375 }
376 const G4MaterialCutsCouple* ghostMCCouple =
378 ->GetMaterialCutsCouple(ghostMaterial,prodCuts);
379 if(ghostMCCouple)
380 {
381 realWorldStepPoint->SetMaterial(ghostMaterial);
382 realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
383 *(fpHyperStep->GetPostStepPoint()) = *(fGhostPostStepPoint);
384 fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
385 fpHyperStep->GetPostStepPoint()->SetMaterialCutsCouple(ghostMCCouple);
386 }
387 else
388 {
389 G4cout << "!!! MaterialCutsCouple is not found for "
390 << ghostMaterial->GetName() << "." << G4endl
391 << " Material in real world ("
392 << realWorldStepPoint->GetMaterial()->GetName()
393 << ") is used." << G4endl;
394 }
395 }
396 }
397}
398
400{
401 G4int pdgCode = partDef->GetPDGEncoding();
402 if(pdgCode==0)
403 {
404 G4String partName = partDef->GetParticleName();
405 if(partName=="opticalphoton") return false;
406 if(partName=="geantino") return false;
407 if(partName=="chargedgeantino") return false;
408 }
409 else
410 {
411 if(pdgCode==22) return false; // gamma
412 if(pdgCode==11) return false; // electron
413 if(pdgCode==2212) return false; // proton
414 if(pdgCode==-12) return false; // anti_nu_e
415 if(pdgCode==12) return false; // nu_e
416 if(pdgCode==-14) return false; // anti_nu_mu
417 if(pdgCode==14) return false; // nu_mu
418 if(pdgCode==-16) return false; // anti_nu_tau
419 if(pdgCode==16) return false; // nu_tau
420 }
421 return true;
422}
423
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
G4ForceCondition
@ StronglyForced
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUndefLimited
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
G4StepStatus
Definition: G4StepStatus.hh:51
@ fGeomBoundary
Definition: G4StepStatus.hh:54
@ fWorldBoundary
Definition: G4StepStatus.hh:52
@ fUndefined
Definition: G4StepStatus.hh:66
@ fPostStepDoItProc
Definition: G4StepStatus.hh:60
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
Definition: G4Material.hh:177
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
static const G4Step * GetHyperStep()
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ParallelWorldProcess(const G4String &processName="ParaWorld", G4ProcessType theType=fParameterisation)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
void SetParallelWorld(G4String parallelWorldName)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
G4bool IsAtRestRequired(G4ParticleDefinition *)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
const G4String & GetParticleName() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ProductionCuts * GetProductionCuts() const
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetMaterial(G4Material *)
G4StepStatus GetStepStatus() const
void SetStepStatus(const G4StepStatus aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
const G4TouchableHandle & GetTouchableHandle() const
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VSensitiveDetector * GetSensitiveDetector() const
Definition: G4Step.hh:78
G4SteppingControl GetControlFlag() const
G4Track * GetTrack() const
void SetStepLength(G4double value)
void SetNonIonizingEnergyDeposit(G4double value)
G4StepPoint * GetPreStepPoint() const
G4double GetNonIonizingEnergyDeposit() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
void SetControlFlag(G4SteppingControl StepControlFlag)
void SetTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const
void SetTrack(G4Track *value)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
const G4Step * GetStep() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
virtual void Initialize(const G4Track &)
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
G4bool Hit(G4Step *aStep)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_MAX
Definition: templates.hh:83